xref: /petsc/src/mat/impls/cufft/cufft.cu (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1e133240eSMatthew G Knepley 
2e133240eSMatthew G Knepley /*
3e133240eSMatthew G Knepley     Provides an interface to the CUFFT package.
4*c4762a1bSJed Brown     Testing examples can be found in ~src/mat/tests
5e133240eSMatthew G Knepley */
6e133240eSMatthew G Knepley 
7af0996ceSBarry Smith #include <petsc/private/matimpl.h>          /*I "petscmat.h" I*/
8e133240eSMatthew G Knepley EXTERN_C_BEGIN
9e133240eSMatthew G Knepley #include <cuda.h>
10e133240eSMatthew G Knepley #include <cuda_runtime.h>
11e133240eSMatthew G Knepley #include <cufft.h>
12e133240eSMatthew G Knepley EXTERN_C_END
13e133240eSMatthew G Knepley 
14e133240eSMatthew G Knepley typedef struct {
15e133240eSMatthew G Knepley   PetscInt     ndim;
16e133240eSMatthew G Knepley   PetscInt     *dim;
17e133240eSMatthew G Knepley   cufftHandle  p_forward, p_backward;
18e133240eSMatthew G Knepley   cufftComplex *devArray;
19e133240eSMatthew G Knepley } Mat_CUFFT;
20e133240eSMatthew G Knepley 
21e133240eSMatthew G Knepley PetscErrorCode MatMult_SeqCUFFT(Mat A, Vec x, Vec y)
22e133240eSMatthew G Knepley {
23e133240eSMatthew G Knepley   Mat_CUFFT      *cufft    = (Mat_CUFFT*) A->data;
24e133240eSMatthew G Knepley   cufftComplex   *devArray = cufft->devArray;
25e133240eSMatthew G Knepley   PetscInt       ndim      = cufft->ndim, *dim = cufft->dim;
26e133240eSMatthew G Knepley   PetscScalar    *x_array, *y_array;
27e133240eSMatthew G Knepley   cufftResult    result;
28e133240eSMatthew G Knepley   PetscErrorCode ierr;
29e133240eSMatthew G Knepley 
30e133240eSMatthew G Knepley   PetscFunctionBegin;
31e133240eSMatthew G Knepley   ierr = VecGetArray(x, &x_array);CHKERRQ(ierr);
32e133240eSMatthew G Knepley   ierr = VecGetArray(y, &y_array);CHKERRQ(ierr);
33e133240eSMatthew G Knepley   if (!cufft->p_forward) {
34e133240eSMatthew G Knepley     cufftResult result;
35e133240eSMatthew G Knepley     /* create a plan, then execute it */
36e133240eSMatthew G Knepley     switch (ndim) {
37e133240eSMatthew G Knepley     case 1:
38e133240eSMatthew G Knepley       result = cufftPlan1d(&cufft->p_forward, dim[0], CUFFT_C2C, 1);CHKERRQ(result != CUFFT_SUCCESS);
39e133240eSMatthew G Knepley       break;
40e133240eSMatthew G Knepley     case 2:
41e133240eSMatthew G Knepley       result = cufftPlan2d(&cufft->p_forward, dim[0], dim[1], CUFFT_C2C);CHKERRQ(result != CUFFT_SUCCESS);
42e133240eSMatthew G Knepley       break;
43e133240eSMatthew G Knepley     case 3:
44e133240eSMatthew G Knepley       result = cufftPlan3d(&cufft->p_forward, dim[0], dim[1], dim[2], CUFFT_C2C);CHKERRQ(result != CUFFT_SUCCESS);
45e133240eSMatthew G Knepley       break;
46e133240eSMatthew G Knepley     default:
47e133240eSMatthew G Knepley       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot create plan for %d-dimensional transform", ndim);
48e133240eSMatthew G Knepley     }
49e133240eSMatthew G Knepley   }
50e133240eSMatthew G Knepley   /* transfer to GPU memory */
51e133240eSMatthew G Knepley   cudaMemcpy(devArray, x_array, sizeof(cufftComplex)*dim[ndim], cudaMemcpyHostToDevice);
52e133240eSMatthew G Knepley   /* execute transform */
53e133240eSMatthew G Knepley   result = cufftExecC2C(cufft->p_forward, devArray, devArray, CUFFT_FORWARD);CHKERRQ(result != CUFFT_SUCCESS);
54e133240eSMatthew G Knepley   /* transfer from GPU memory */
55e133240eSMatthew G Knepley   cudaMemcpy(y_array, devArray, sizeof(cufftComplex)*dim[ndim], cudaMemcpyDeviceToHost);
56e133240eSMatthew G Knepley   ierr = VecRestoreArray(y, &y_array);CHKERRQ(ierr);
57e133240eSMatthew G Knepley   ierr = VecRestoreArray(x, &x_array);CHKERRQ(ierr);
58e133240eSMatthew G Knepley   PetscFunctionReturn(0);
59e133240eSMatthew G Knepley }
60e133240eSMatthew G Knepley 
61e133240eSMatthew G Knepley PetscErrorCode MatMultTranspose_SeqCUFFT(Mat A, Vec x, Vec y)
62e133240eSMatthew G Knepley {
63e133240eSMatthew G Knepley   Mat_CUFFT      *cufft    = (Mat_CUFFT*) A->data;
64e133240eSMatthew G Knepley   cufftComplex   *devArray = cufft->devArray;
65e133240eSMatthew G Knepley   PetscInt       ndim      = cufft->ndim, *dim = cufft->dim;
66e133240eSMatthew G Knepley   PetscScalar    *x_array, *y_array;
67e133240eSMatthew G Knepley   cufftResult    result;
68e133240eSMatthew G Knepley   PetscErrorCode ierr;
69e133240eSMatthew G Knepley 
70e133240eSMatthew G Knepley   PetscFunctionBegin;
71e133240eSMatthew G Knepley   ierr = VecGetArray(x, &x_array);CHKERRQ(ierr);
72e133240eSMatthew G Knepley   ierr = VecGetArray(y, &y_array);CHKERRQ(ierr);
73e133240eSMatthew G Knepley   if (!cufft->p_backward) {
74e133240eSMatthew G Knepley     /* create a plan, then execute it */
75e133240eSMatthew G Knepley     switch (ndim) {
76e133240eSMatthew G Knepley     case 1:
77e133240eSMatthew G Knepley       result = cufftPlan1d(&cufft->p_backward, dim[0], CUFFT_C2C, 1);CHKERRQ(result != CUFFT_SUCCESS);
78e133240eSMatthew G Knepley       break;
79e133240eSMatthew G Knepley     case 2:
80e133240eSMatthew G Knepley       result = cufftPlan2d(&cufft->p_backward, dim[0], dim[1], CUFFT_C2C);CHKERRQ(result != CUFFT_SUCCESS);
81e133240eSMatthew G Knepley       break;
82e133240eSMatthew G Knepley     case 3:
83e133240eSMatthew G Knepley       result = cufftPlan3d(&cufft->p_backward, dim[0], dim[1], dim[2], CUFFT_C2C);CHKERRQ(result != CUFFT_SUCCESS);
84e133240eSMatthew G Knepley       break;
85e133240eSMatthew G Knepley     default:
86e133240eSMatthew G Knepley       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot create plan for %d-dimensional transform", ndim);
87e133240eSMatthew G Knepley     }
88e133240eSMatthew G Knepley   }
89e133240eSMatthew G Knepley   /* transfer to GPU memory */
90e133240eSMatthew G Knepley   cudaMemcpy(devArray, x_array, sizeof(cufftComplex)*dim[ndim], cudaMemcpyHostToDevice);
91e133240eSMatthew G Knepley   /* execute transform */
92e133240eSMatthew G Knepley   result = cufftExecC2C(cufft->p_forward, devArray, devArray, CUFFT_INVERSE);CHKERRQ(result != CUFFT_SUCCESS);
93e133240eSMatthew G Knepley   /* transfer from GPU memory */
94e133240eSMatthew G Knepley   cudaMemcpy(y_array, devArray, sizeof(cufftComplex)*dim[ndim], cudaMemcpyDeviceToHost);
95e133240eSMatthew G Knepley   ierr = VecRestoreArray(y, &y_array);CHKERRQ(ierr);
96e133240eSMatthew G Knepley   ierr = VecRestoreArray(x, &x_array);CHKERRQ(ierr);
97e133240eSMatthew G Knepley   PetscFunctionReturn(0);
98e133240eSMatthew G Knepley }
99e133240eSMatthew G Knepley 
100e133240eSMatthew G Knepley PetscErrorCode MatDestroy_SeqCUFFT(Mat A)
101e133240eSMatthew G Knepley {
102e133240eSMatthew G Knepley   Mat_CUFFT      *cufft = (Mat_CUFFT*) A->data;
103e133240eSMatthew G Knepley   cufftResult    result;
104e133240eSMatthew G Knepley   PetscErrorCode ierr;
105e133240eSMatthew G Knepley 
106e133240eSMatthew G Knepley   PetscFunctionBegin;
107e133240eSMatthew G Knepley   ierr = PetscFree(cufft->dim);CHKERRQ(ierr);
108e133240eSMatthew G Knepley   if (cufft->p_forward)  {result = cufftDestroy(cufft->p_forward);CHKERRQ(result != CUFFT_SUCCESS);}
109e133240eSMatthew G Knepley   if (cufft->p_backward) {result = cufftDestroy(cufft->p_backward);CHKERRQ(result != CUFFT_SUCCESS);}
110e133240eSMatthew G Knepley   cudaFree(cufft->devArray);
111bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
112bf0cc555SLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
113e133240eSMatthew G Knepley   PetscFunctionReturn(0);
114e133240eSMatthew G Knepley }
115e133240eSMatthew G Knepley 
116e133240eSMatthew G Knepley /*@
117e133240eSMatthew G Knepley   MatCreateSeqCUFFT - Creates a matrix object that provides sequential FFT via the external package CUFFT
118e133240eSMatthew G Knepley 
119d083f849SBarry Smith   Collective
120e133240eSMatthew G Knepley 
121e133240eSMatthew G Knepley   Input Parameters:
122e133240eSMatthew G Knepley + comm - MPI communicator, set to PETSC_COMM_SELF
123e133240eSMatthew G Knepley . ndim - the ndim-dimensional transform
124e133240eSMatthew G Knepley - dim  - array of size ndim, dim[i] contains the vector length in the i-dimension
125e133240eSMatthew G Knepley 
126e133240eSMatthew G Knepley   Output Parameter:
127e133240eSMatthew G Knepley . A - the matrix
128e133240eSMatthew G Knepley 
129e133240eSMatthew G Knepley   Options Database Keys:
130e133240eSMatthew G Knepley . -mat_cufft_plannerflags - set CUFFT planner flags
131e133240eSMatthew G Knepley 
132e133240eSMatthew G Knepley   Level: intermediate
133e133240eSMatthew G Knepley @*/
1347087cfbeSBarry Smith PetscErrorCode  MatCreateSeqCUFFT(MPI_Comm comm, PetscInt ndim, const PetscInt dim[], Mat *A)
135e133240eSMatthew G Knepley {
136e133240eSMatthew G Knepley   Mat_CUFFT      *cufft;
137e133240eSMatthew G Knepley   PetscInt       m, d;
138e133240eSMatthew G Knepley   PetscErrorCode ierr;
139e133240eSMatthew G Knepley 
140e133240eSMatthew G Knepley   PetscFunctionBegin;
141e133240eSMatthew G Knepley   if (ndim < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_USER, "ndim %d must be > 0", ndim);
142e133240eSMatthew G Knepley   ierr = MatCreate(comm, A);CHKERRQ(ierr);
143e133240eSMatthew G Knepley   m    = 1;
144e133240eSMatthew G Knepley   for (d = 0; d < ndim; ++d) {
145e133240eSMatthew G Knepley     if (dim[d] < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_USER, "dim[%d]=%d must be > 0", d, dim[d]);
146e133240eSMatthew G Knepley     m *= dim[d];
147e133240eSMatthew G Knepley   }
148e133240eSMatthew G Knepley   ierr = MatSetSizes(*A, m, m, m, m);CHKERRQ(ierr);
149e133240eSMatthew G Knepley   ierr = PetscObjectChangeTypeName((PetscObject)*A, MATSEQCUFFT);CHKERRQ(ierr);
150e133240eSMatthew G Knepley 
151b00a9115SJed Brown   ierr       = PetscNewLog(*A,&cufft);CHKERRQ(ierr);
152e133240eSMatthew G Knepley   (*A)->data = (void*) cufft;
153854ce69bSBarry Smith   ierr       = PetscMalloc1(ndim+1, &cufft->dim);CHKERRQ(ierr);
154580bdb30SBarry Smith   ierr       = PetscArraycpy(cufft->dim, dim, ndim);CHKERRQ(ierr);
1552205254eSKarl Rupp 
156e133240eSMatthew G Knepley   cufft->ndim       = ndim;
157e133240eSMatthew G Knepley   cufft->p_forward  = 0;
158e133240eSMatthew G Knepley   cufft->p_backward = 0;
159e133240eSMatthew G Knepley   cufft->dim[ndim]  = m;
160e133240eSMatthew G Knepley 
161e133240eSMatthew G Knepley   /* GPU memory allocation */
162e133240eSMatthew G Knepley   cudaMalloc((void**) &cufft->devArray, sizeof(cufftComplex)*m);
163e133240eSMatthew G Knepley 
164e133240eSMatthew G Knepley   (*A)->ops->mult          = MatMult_SeqCUFFT;
165e133240eSMatthew G Knepley   (*A)->ops->multtranspose = MatMultTranspose_SeqCUFFT;
166e133240eSMatthew G Knepley   (*A)->assembled          = PETSC_TRUE;
167e133240eSMatthew G Knepley   (*A)->ops->destroy       = MatDestroy_SeqCUFFT;
168e133240eSMatthew G Knepley 
169e133240eSMatthew G Knepley   /* get runtime options */
170e133240eSMatthew G Knepley   ierr = PetscOptionsBegin(comm, ((PetscObject)(*A))->prefix, "CUFFT Options", "Mat");CHKERRQ(ierr);
171e133240eSMatthew G Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
172e133240eSMatthew G Knepley   PetscFunctionReturn(0);
173e133240eSMatthew G Knepley }
174