hipFFT API Usage#

This section describes usage of the hipFFT library API. The hipFFT API follows the NVIDIA cuFFT API.

Types#

There are a few data structures that are internal to the library. The pointer types to these structures are given below. The user would need to use these types to create handles and pass them between different library functions.

HIPFFT_FORWARD#

Perform a forward FFT.

HIPFFT_BACKWARD#

Perform a backward/inverse FFT.

Warning

doxygenenum: Cannot find enum “hipfftType_t” in doxygen xml output for project “hipFFT 1.0.17 Documentation” from directory: /home/docs/checkouts/readthedocs.org/user_builds/advanced-micro-devices-hipfft/checkouts/latest/docs/doxygen/xml

typedef struct hipfftHandle_t *hipfftHandle#

Warning

doxygenenum: Cannot find enum “hipfftResult_t” in doxygen xml output for project “hipFFT 1.0.17 Documentation” from directory: /home/docs/checkouts/readthedocs.org/user_builds/advanced-micro-devices-hipfft/checkouts/latest/docs/doxygen/xml

Simple plans#

These planning routines allocate a plan for you. If execution of the plan requires a work buffer, it will be created (and destroyed) automatically.

hipfftResult hipfftPlan1d(hipfftHandle *plan, int nx, hipfftType type, int batch)#

Create a new one-dimensional FFT plan.

Allocate and initialize a new one-dimensional FFT plan.

Parameters:
  • plan[out] Pointer to the FFT plan handle.

  • nx[in] FFT length.

  • type[in] FFT type.

  • batch[in] Number of batched transforms to compute.

hipfftResult hipfftPlan2d(hipfftHandle *plan, int nx, int ny, hipfftType type)#

Create a new two-dimensional FFT plan.

Allocate and initialize a new two-dimensional FFT plan. Two-dimensional data should be stored in C ordering (row-major format), so that indexes in y-direction (j index) vary the fastest.

Parameters:
  • plan[out] Pointer to the FFT plan handle.

  • nx[in] Number of elements in the x-direction (slow index).

  • ny[in] Number of elements in the y-direction (fast index).

  • type[in] FFT type.

hipfftResult hipfftPlan3d(hipfftHandle *plan, int nx, int ny, int nz, hipfftType type)#

Create a new three-dimensional FFT plan.

Allocate and initialize a new three-dimensional FFT plan. Three-dimensional data should be stored in C ordering (row-major format), so that indexes in z-direction (k index) vary the fastest.

Parameters:
  • plan[out] Pointer to the FFT plan handle.

  • nx[in] Number of elements in the x-direction (slowest index).

  • ny[in] Number of elements in the y-direction.

  • nz[in] Number of elements in the z-direction (fastest index).

  • type[in] FFT type.

User managed simple plans#

These planning routines assume that you have allocated a plan (hipfftHandle) yourself; and that you will manage a work area as well.

If you want to manage your own work buffer… XXX

hipfftResult hipfftCreate(hipfftHandle *plan)#

Allocate a new plan.

hipfftResult hipfftDestroy(hipfftHandle plan)#

Destroy and deallocate an existing plan.

hipfftResult hipfftSetAutoAllocation(hipfftHandle plan, int autoAllocate)#

Set the plan’s auto-allocation flag. The plan will allocate its own workarea.

Parameters:
  • plan[in] Pointer to the FFT plan.

  • autoAllocate[in] 0 to disable auto-allocation, non-zero to enable.

hipfftResult hipfftMakePlan1d(hipfftHandle plan, int nx, hipfftType type, int batch, size_t *workSize)#

Initialize a new one-dimensional FFT plan.

Assumes that the plan has been created already, and modifies the plan associated with the plan handle.

Parameters:
  • plan[in] Handle of the FFT plan.

  • nx[in] FFT length.

  • type[in] FFT type.

  • batch[in] Number of batched transforms to compute.

  • workSize[out] Pointer to work area size (returned value).

hipfftResult hipfftMakePlan2d(hipfftHandle plan, int nx, int ny, hipfftType type, size_t *workSize)#

Initialize a new two-dimensional FFT plan.

Assumes that the plan has been created already, and modifies the plan associated with the plan handle. Two-dimensional data should be stored in C ordering (row-major format), so that indexes in y-direction (j index) vary the fastest.

Parameters:
  • plan[in] Handle of the FFT plan.

  • nx[in] Number of elements in the x-direction (slow index).

  • ny[in] Number of elements in the y-direction (fast index).

  • type[in] FFT type.

  • workSize[out] Pointer to work area size (returned value).

hipfftResult hipfftMakePlan3d(hipfftHandle plan, int nx, int ny, int nz, hipfftType type, size_t *workSize)#

Initialize a new two-dimensional FFT plan.

Assumes that the plan has been created already, and modifies the plan associated with the plan handle. Three-dimensional data should be stored in C ordering (row-major format), so that indexes in z-direction (k index) vary the fastest.

Parameters:
  • plan[in] Handle of the FFT plan.

  • nx[in] Number of elements in the x-direction (slowest index).

  • ny[in] Number of elements in the y-direction.

  • nz[in] Number of elements in the z-direction (fastest index).

  • type[in] FFT type.

  • workSize[out] Pointer to work area size (returned value).

Advanced plans#

hipfftResult hipfftMakePlanMany(hipfftHandle plan, int rank, int *n, int *inembed, int istride, int idist, int *onembed, int ostride, int odist, hipfftType type, int batch, size_t *workSize)#

Initialize a new batched rank-dimensional FFT plan with advanced data layout.

Assumes that the plan has been created already, and modifies the plan associated with the plan handle. The number of elements to transform in each direction of the input data in the FFT plan is specified in n.

The batch parameter tells hipFFT how many transforms to perform. The distance between the first elements of two consecutive batches of the input and output data are specified with the idist and odist parameters.

The inembed and onembed parameters define the input and output data layouts. The number of elements in the data is assumed to be larger than the number of elements in the transform. Strided data layouts are also supported. Strides along the fastest direction in the input and output data are specified via the istride and ostride parameters.

If both inembed and onembed parameters are set to NULL, all the advanced data layout parameters are ignored and reverted to default values, i.e., the batched transform is performed with non-strided data access and the number of data/transform elements are assumed to be

equivalent.

Parameters:
  • plan[out] Pointer to the FFT plan handle.

  • rank[in] Dimension of transform (1, 2, or 3).

  • n[in] Number of elements to transform in the x/y/z directions.

  • inembed[in] Number of elements in the input data in the x/y/z directions.

  • istride[in] Distance between two successive elements in the input data.

  • idist[in] Distance between input batches.

  • onembed[in] Number of elements in the output data in the x/y/z directions.

  • ostride[in] Distance between two successive elements in the output data.

  • odist[in] Distance between output batches.

  • type[in] FFT type.

  • batch[in] Number of batched transforms to perform.

  • workSize[out] Pointer to work area size (returned value).

hipfftResult hipfftXtMakePlanMany(hipfftHandle plan, int rank, long long int *n, long long int *inembed, long long int istride, long long int idist, hipDataType inputType, long long int *onembed, long long int ostride, long long int odist, hipDataType outputType, long long int batch, size_t *workSize, hipDataType executionType)#

Initialize a batched rank-dimensional FFT plan with advanced data layout and specified input, output, execution data types.

Assumes that the plan has been created already, and modifies the plan associated with the plan handle. The number of elements to transform in each direction of the input data in the FFT plan is specified in n.

The batch parameter tells hipFFT how many transforms to perform. The distance between the first elements of two consecutive batches of the input and output data are specified with the idist and odist parameters.

The inembed and onembed parameters define the input and output data layouts. The number of elements in the data is assumed to be larger than the number of elements in the transform. Strided data layouts are also supported. Strides along the fastest direction in the input and output data are specified via the istride and ostride parameters.

If both inembed and onembed parameters are set to NULL, all the advanced data layout parameters are ignored and reverted to default values, i.e., the batched transform is performed with non-strided data access and the number of data/transform elements are assumed to be

equivalent.

The inputType, outputType, executionType parameters specify the data types (precision, and whether the data is real-valued or complex-valued) of the transform input, output, and internal representation respectively. Currently, the precision of all three parameters must match, and the execution type must always be complex-valued. At least one of inputType and outputType must be complex. A half-precision transform can be requested by using either the HIP_R_16F or HIP_C_16F types.

Parameters:
  • plan[out] Pointer to the FFT plan handle.

  • rank[in] Dimension of transform (1, 2, or 3).

  • n[in] Number of elements to transform in the x/y/z directions.

  • inembed[in] Number of elements in the input data in the x/y/z directions.

  • istride[in] Distance between two successive elements in the input data.

  • idist[in] Distance between input batches.

  • inputType[in] Format of FFT input.

  • onembed[in] Number of elements in the output data in the x/y/z directions.

  • ostride[in] Distance between two successive elements in the output data.

  • odist[in] Distance between output batches.

  • outputType[in] Format of FFT output.

  • batch[in] Number of batched transforms to perform.

  • workSize[out] Pointer to work area size (returned value).

  • executionType[in] Internal data format used by the library during computation.

Estimating work area sizes#

These call return estimates of the work area required to support a plan generated with the same parameters (either with the simple or extensible API). Callers who choose to manage work area allocation within their application must use this call after plan generation, and after any hipfftSet*() calls subsequent to plan generation, if those calls might alter the required work space size.

hipfftResult hipfftEstimate1d(int nx, hipfftType type, int batch, size_t *workSize)#

Return an estimate of the work area size required for a 1D plan.

Parameters:
  • nx[in] Number of elements in the x-direction.

  • type[in] FFT type.

  • batch[in] Number of batched transforms to perform.

  • workSize[out] Pointer to work area size (returned value).

hipfftResult hipfftEstimate2d(int nx, int ny, hipfftType type, size_t *workSize)#

Return an estimate of the work area size required for a 2D plan.

Parameters:
  • nx[in] Number of elements in the x-direction.

  • ny[in] Number of elements in the y-direction.

  • type[in] FFT type.

  • workSize[out] Pointer to work area size (returned value).

hipfftResult hipfftEstimate3d(int nx, int ny, int nz, hipfftType type, size_t *workSize)#

Return an estimate of the work area size required for a 3D plan.

Parameters:
  • nx[in] Number of elements in the x-direction.

  • ny[in] Number of elements in the y-direction.

  • nz[in] Number of elements in the z-direction.

  • type[in] FFT type.

  • workSize[out] Pointer to work area size (returned value).

hipfftResult hipfftEstimateMany(int rank, int *n, int *inembed, int istride, int idist, int *onembed, int ostride, int odist, hipfftType type, int batch, size_t *workSize)#

Return an estimate of the work area size required for a rank-dimensional plan.

Parameters:
  • rank[in] Dimension of FFT transform (1, 2, or 3).

  • n[in] Number of elements in the x/y/z directions.

  • inembed[in]

  • istride[in]

  • idist[in] Distance between input batches.

  • onembed[in]

  • ostride[in]

  • odist[in] Distance between output batches.

  • type[in] FFT type.

  • batch[in] Number of batched transforms to perform.

  • workSize[out] Pointer to work area size (returned value).

Accurate work area sizes#

After plan generation is complete, an accurate work area size can be obtained with these routines.

hipfftResult hipfftGetSize1d(hipfftHandle plan, int nx, hipfftType type, int batch, size_t *workSize)#

Return size of the work area size required for a 1D plan.

Parameters:
  • plan[in] Pointer to the FFT plan.

  • nx[in] Number of elements in the x-direction.

  • type[in] FFT type.

  • batch[in] Number of batched transforms to perform.

  • workSize[out] Pointer to work area size (returned value).

hipfftResult hipfftGetSize2d(hipfftHandle plan, int nx, int ny, hipfftType type, size_t *workSize)#

Return size of the work area size required for a 2D plan.

Parameters:
  • plan[in] Pointer to the FFT plan.

  • nx[in] Number of elements in the x-direction.

  • ny[in] Number of elements in the y-direction.

  • type[in] FFT type.

  • workSize[out] Pointer to work area size (returned value).

hipfftResult hipfftGetSize3d(hipfftHandle plan, int nx, int ny, int nz, hipfftType type, size_t *workSize)#

Return size of the work area size required for a 3D plan.

Parameters:
  • plan[in] Pointer to the FFT plan.

  • nx[in] Number of elements in the x-direction.

  • ny[in] Number of elements in the y-direction.

  • nz[in] Number of elements in the z-direction.

  • type[in] FFT type.

  • workSize[out] Pointer to work area size (returned value).

hipfftResult hipfftGetSizeMany(hipfftHandle plan, int rank, int *n, int *inembed, int istride, int idist, int *onembed, int ostride, int odist, hipfftType type, int batch, size_t *workSize)#

Return size of the work area size required for a rank-dimensional plan.

Parameters:
  • plan[in] Pointer to the FFT plan.

  • rank[in] Dimension of FFT transform (1, 2, or 3).

  • n[in] Number of elements in the x/y/z directions.

  • inembed[in]

  • istride[in]

  • idist[in] Distance between input batches.

  • onembed[in]

  • ostride[in]

  • odist[in] Distance between output batches.

  • type[in] FFT type.

  • batch[in] Number of batched transforms to perform.

  • workSize[out] Pointer to work area size (returned value).

hipfftResult hipfftXtGetSizeMany(hipfftHandle plan, int rank, long long int *n, long long int *inembed, long long int istride, long long int idist, hipDataType inputType, long long int *onembed, long long int ostride, long long int odist, hipDataType outputType, long long int batch, size_t *workSize, hipDataType executionType)#

Return size of the work area size required for a rank-dimensional plan, with specified input, output, execution data types.

See hipfftXtMakePlanMany for restrictions on inputType, outputType, executionType parameters.

Parameters:
  • plan[in] Pointer to the FFT plan.

  • rank[in] Dimension of FFT transform (1, 2, or 3).

  • n[in] Number of elements in the x/y/z directions.

  • inembed[in] Number of elements in the input data in the x/y/z directions.

  • istride[in] Distance between two successive elements in the input data.

  • idist[in] Distance between input batches.

  • inputType[in] Format of FFT input.

  • onembed[in] Number of elements in the output data in the x/y/z directions.

  • ostride[in] Distance between two successive elements in the output data.

  • odist[in] Distance between output batches.

  • outputType[in] Format of FFT output.

  • batch[in] Number of batched transforms to perform.

  • workSize[out] Pointer to work area size (returned value).

  • executionType[in] Internal data format used by the library during computation.

Executing plans#

Once you have created an FFT plan, you can execute it using one of the hipfftExec* functions.

For real-to-complex transforms, the output buffer XXX

For complex-to-real transforms, the output buffer XXX

hipfftResult hipfftExecC2C(hipfftHandle plan, hipfftComplex *idata, hipfftComplex *odata, int direction)#

Execute a (float) complex-to-complex FFT.

If the input and output buffers are equal, an in-place transform is performed.

Parameters:
  • plan – The FFT plan.

  • idata – Input data (on device).

  • odata – Output data (on device).

  • direction – Either HIPFFT_FORWARD or HIPFFT_BACKWARD.

hipfftResult hipfftExecR2C(hipfftHandle plan, hipfftReal *idata, hipfftComplex *odata)#

Execute a (float) real-to-complex FFT.

If the input and output buffers are equal, an in-place transform is performed.

Parameters:
  • plan – The FFT plan.

  • idata – Input data (on device).

  • odata – Output data (on device).

hipfftResult hipfftExecC2R(hipfftHandle plan, hipfftComplex *idata, hipfftReal *odata)#

Execute a (float) complex-to-real FFT.

If the input and output buffers are equal, an in-place transform is performed.

Parameters:
  • plan – The FFT plan.

  • idata – Input data (on device).

  • odata – Output data (on device).

hipfftResult hipfftExecZ2Z(hipfftHandle plan, hipfftDoubleComplex *idata, hipfftDoubleComplex *odata, int direction)#

Execute a (double) complex-to-complex FFT.

If the input and output buffers are equal, an in-place transform is performed.

Parameters:
  • plan – The FFT plan.

  • idata – Input data (on device).

  • odata – Output data (on device).

  • direction – Either HIPFFT_FORWARD or HIPFFT_BACKWARD.

hipfftResult hipfftExecD2Z(hipfftHandle plan, hipfftDoubleReal *idata, hipfftDoubleComplex *odata)#

Execute a (double) real-to-complex FFT.

If the input and output buffers are equal, an in-place transform is performed.

Parameters:
  • plan – The FFT plan.

  • idata – Input data (on device).

  • odata – Output data (on device).

hipfftResult hipfftExecZ2D(hipfftHandle plan, hipfftDoubleComplex *idata, hipfftDoubleReal *odata)#

Execute a (double) complex-to-real FFT.

If the input and output buffers are equal, an in-place transform is performed.

Parameters:
  • plan – The FFT plan.

  • idata – Input data (on device).

  • odata – Output data (on device).

hipfftResult hipfftXtExec(hipfftHandle plan, void *input, void *output, int direction)#

Execute an FFT plan for any precision and type.

An in-place transform is performed if the input and output pointers have the same value.

The direction parameter is ignored if for real-to-complex and complex-to-real transforms, as the direction is already implied by the data types.

Parameters:
  • plan[in] Pointer to the FFT plan.

  • input[in] Pointer to input data for the transform.

  • output[in] Pointer to output data for the transform.

  • direction[in] Either HIPFFT_FORWARD or HIPFFT_BACKWARD.

HIP graph support for hipFFT#

hipFFT supports capturing kernels launched during FFT execution into HIP graph nodes. This way, users can capture FFT execution, along with other work, into a HIP graph and launch the work in the graph multiple times.

The following hipFFT APIs can be used with graph capture:

Note that each launch of a HIP graph will provide the same arguments to the kernels in the graph. In particular, this implies that all of the parameters to the above APIs remain valid while the HIP graph is in use:

  • The hipFFT plan

  • The input and output buffers

hipFFT does not support capturing work performed by other API functions aside from those listed above.

Callbacks#

hipfftResult hipfftXtSetCallback(hipfftHandle plan, void **callbacks, hipfftXtCallbackType cbtype, void **callbackData)#

Set a callback on a plan.

Set either a load or store callback to run with a plan. The type of callback is specified with the ‘cbtype’ parameter. An array ofcallback and callback data pointers must be given - one per device executing the plan.

Parameters:
  • plan[in] The FFT plan.

  • callbacks[in] Array of callback function pointers.

  • cbtype[in] Type of callback being set.

  • callbackData[in] Array of callback function data pointers

hipfftResult hipfftXtClearCallback(hipfftHandle plan, hipfftXtCallbackType cbtype)#

Remove a callback from a plan.

Remove a previously-set callback from a plan.

Parameters:
  • plan[in] The FFT plan.

  • cbtype[in] Type of callback being removed.

hipfftResult hipfftXtSetCallbackSharedSize(hipfftHandle plan, hipfftXtCallbackType cbtype, size_t sharedSize)#

Set shared memory size for callback.

Set shared memory required for a callback. The callback of the specified type must have already been set on the plan.

Parameters:
  • plan[in] The FFT plan.

  • cbtype[in] Type of callback being modified.

  • sharedSize[in] Amount of shared memory required, in bytes.

Single-process Multi-GPU Transforms#

hipFFT offers experimental single-process multi-GPU transforms.

Multiple devices can be associated to a hipfftHandle using hipfftXtSetGPUs(). Once a plan is associated to multiple GPUs, hipLibXtDesc is used to pass multiple GPU buffers to the plan for execution. These buffers are allocated via hipfftXtMalloc(), and free’d with hipfftXtFree(). The function hipfftXtMemcpy() allows one to move data to or from a hipLibXtDesc and a contiguous host buffer, or between two hipLibXtDesc s.

Execution is performed with the appropriate hipfftXtExecDescriptor()

hipfftResult hipfftXtSetGPUs(hipfftHandle plan, int count, int *gpus)#

Set multiple GPUs on a plan.

Instructs hipFFT to use multiple GPUs for a plan.

This function must be called after the plan is allocated using hipfftCreate, but before the plan is initialized by any of the “MakePlan” functions. Therefore, API functions that combine creation and initialization (hipfftPlan1d, hipfftPlan2d, hipfftPlan3d, and hipfftPlanMany) cannot use multiple GPUs.

Warning

Experimental

Parameters:
  • plan[inout]

  • count[in] length gpus array

  • gpus[in] array of ints specifying deviceIDs

struct hipXtDesc_t#

Struct for single-process multi-GPU transform.

This struct holds pointers to device memory, including the device the memory resides on and the size of each block of memory.

Warning

Experimental

struct hipLibXtDesc_t#

Struct for single-process multi-GPU transform.

This struct holds hipXtDesc_t structures that define blocks of memory for use in a transform.

Warning

Experimental

hipfftResult hipfftXtMalloc(hipfftHandle plan, hipLibXtDesc **desc, hipfftXtSubFormat format)#

Allocate memory on multiple devices.

Allocate memory on multiple devices for the specified plan. Returns a hipLibXtDesc_t descriptor which includes pointers to the allocated memory, devices that memory resides on, and sizes allocated.

The subformat indicates whether the memory will be used for FFT input or output.

The memory must be freed by calling hipfftXtFree.

Warning

Experimental

hipfftResult hipfftXtFree(hipLibXtDesc *desc)#

Free memory allocated by hipfftXtMalloc.

Warning

Experimental

hipfftResult hipfftXtMemcpy(hipfftHandle plan, void *dest, void *src, hipfftXtCopyType type)#

Copy data to/from hipLibXtDesc_t descriptors.

Copy data according to the hipfftXtCopyType_t

  • HIPFFT_COPY_HOST_TO_DEVICE: dest points to a hipLibXtDesc_t structure that describes multi-device memory layout. src points to a host memory buffer.

  • HIPFFT_COPY_DEVICE_TO_HOST: src points to a hipLibXtDesc_t structure that describes multi-device memory layout. dest points to a host memory buffer.

  • HIPFFT_COPY_DEVICE_TO_DEVICE: Both dest and src point to a hipLibXtDesc_t structure that describes multi-device memory layout. The two structures must describe memory with the same number of devices and memory sizes.

Warning

Experimental

group hipfftXtExecDescriptor

Execute FFTs using hipLibXtDesc_t descriptors. Inputs and outputs are pointers to hipLibXtDesc_t descriptors. In-place transforms are performed by passing the same pointer for input and output.

Warning

Experimental

Functions

hipfftResult hipfftXtExecDescriptorC2C(hipfftHandle plan, hipLibXtDesc *input, hipLibXtDesc *output, int direction)#

Execute single-precision complex-to-complex transform.

hipfftResult hipfftXtExecDescriptorR2C(hipfftHandle plan, hipLibXtDesc *input, hipLibXtDesc *output)#

Execute single-precision real forward transform.

hipfftResult hipfftXtExecDescriptorC2R(hipfftHandle plan, hipLibXtDesc *input, hipLibXtDesc *output)#

Execute single-precision real backward transform.

hipfftResult hipfftXtExecDescriptorZ2Z(hipfftHandle plan, hipLibXtDesc *input, hipLibXtDesc *output, int direction)#

Execute double-precision complex-to-complex transform.

hipfftResult hipfftXtExecDescriptorD2Z(hipfftHandle plan, hipLibXtDesc *input, hipLibXtDesc *output)#

Execute double-precision real forward transform.

hipfftResult hipfftXtExecDescriptorZ2D(hipfftHandle plan, hipLibXtDesc *input, hipLibXtDesc *output)#

Execute double-precision real backward transform.

hipfftResult hipfftXtExecDescriptor(hipfftHandle plan, hipLibXtDesc *input, hipLibXtDesc *output, int direction)#

Execute general transform - types of the descriptors must match the plan.