xref: /petsc/src/mat/impls/fft/fft.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
1dedccee8SHong Zhang /*
2dedccee8SHong Zhang     Provides an interface to the FFT packages.
3dedccee8SHong Zhang */
4dedccee8SHong Zhang 
5c6db04a5SJed Brown #include <../src/mat/impls/fft/fft.h> /*I "petscmat.h" I*/
6dedccee8SHong Zhang 
79371c9d4SSatish Balay PetscErrorCode MatDestroy_FFT(Mat A) {
8dedccee8SHong Zhang   Mat_FFT *fft = (Mat_FFT *)A->data;
9dedccee8SHong Zhang 
10dedccee8SHong Zhang   PetscFunctionBegin;
119566063dSJacob Faibussowitsch   if (fft->matdestroy) PetscCall((fft->matdestroy)(A));
129566063dSJacob Faibussowitsch   PetscCall(PetscFree(fft->dim));
139566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
149566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, NULL));
15dedccee8SHong Zhang   PetscFunctionReturn(0);
16dedccee8SHong Zhang }
17dedccee8SHong Zhang 
180a9977b2SMatthew G. Knepley /*@C
19dedccee8SHong Zhang       MatCreateFFT - Creates a matrix object that provides FFT via an external package
20dedccee8SHong Zhang 
21d083f849SBarry Smith    Collective
22dedccee8SHong Zhang 
23d8d19677SJose E. Roman    Input Parameters:
24dedccee8SHong Zhang +   comm - MPI communicator
25dedccee8SHong Zhang .   ndim - the ndim-dimensional transform
26dedccee8SHong Zhang .   dim - array of size ndim, dim[i] contains the vector length in the i-dimension
2711a5261eSBarry Smith -   type - package type, e.g., `MATFFTW` or `MATSEQCUFFT`
28dedccee8SHong Zhang 
29dedccee8SHong Zhang    Output Parameter:
30dedccee8SHong Zhang .   A  - the matrix
31dedccee8SHong Zhang 
32dedccee8SHong Zhang    Options Database Keys:
3375f45d78SBarry Smith .   -mat_fft_type - set FFT type fft or seqcufft
3475f45d78SBarry Smith 
3511a5261eSBarry Smith    Note: this serves as a base class for all FFT marix classes, currently `MATFFTW` or `MATSEQCUFFT`
36dedccee8SHong Zhang 
37dedccee8SHong Zhang    Level: intermediate
38dedccee8SHong Zhang 
3911a5261eSBarry Smith .seealso: `MATFFTW`, `MATSEQCUFFT`, `MatCreateVecsFFTW()`
40dedccee8SHong Zhang @*/
419371c9d4SSatish Balay PetscErrorCode MatCreateFFT(MPI_Comm comm, PetscInt ndim, const PetscInt dim[], MatType mattype, Mat *A) {
42dedccee8SHong Zhang   PetscMPIInt size;
43dedccee8SHong Zhang   Mat         FFT;
44dedccee8SHong Zhang   PetscInt    N, i;
45dedccee8SHong Zhang   Mat_FFT    *fft;
46dedccee8SHong Zhang 
47dedccee8SHong Zhang   PetscFunctionBegin;
48c0aa6a63SJacob Faibussowitsch   PetscValidIntPointer(dim, 3);
49c0aa6a63SJacob Faibussowitsch   PetscValidPointer(A, 5);
505f80ce2aSJacob Faibussowitsch   PetscCheck(ndim >= 1, comm, PETSC_ERR_USER, "ndim %" PetscInt_FMT " must be > 0", ndim);
519566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
52dedccee8SHong Zhang 
539566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &FFT));
54*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&fft));
55dedccee8SHong Zhang   FFT->data = (void *)fft;
56dedccee8SHong Zhang   N         = 1;
57dedccee8SHong Zhang   for (i = 0; i < ndim; i++) {
585f80ce2aSJacob Faibussowitsch     PetscCheck(dim[i] >= 1, PETSC_COMM_SELF, PETSC_ERR_USER, "dim[%" PetscInt_FMT "]=%" PetscInt_FMT " must be > 0", i, dim[i]);
59dedccee8SHong Zhang     N *= dim[i];
60dedccee8SHong Zhang   }
61dedccee8SHong Zhang 
629566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ndim, &fft->dim));
639566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(fft->dim, dim, ndim));
64dedccee8SHong Zhang 
65dedccee8SHong Zhang   fft->ndim = ndim;
66dedccee8SHong Zhang   fft->n    = PETSC_DECIDE;
67dedccee8SHong Zhang   fft->N    = N;
680298fd71SBarry Smith   fft->data = NULL;
69dedccee8SHong Zhang 
709566063dSJacob Faibussowitsch   PetscCall(MatSetType(FFT, mattype));
7126fbe8dcSKarl Rupp 
72dedccee8SHong Zhang   FFT->ops->destroy = MatDestroy_FFT;
73dedccee8SHong Zhang 
745f80ce2aSJacob Faibussowitsch   /* get runtime options... what options? */
75d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)FFT);
76d0609cedSBarry Smith   PetscOptionsEnd();
77dedccee8SHong Zhang 
78dedccee8SHong Zhang   *A = FFT;
79dedccee8SHong Zhang   PetscFunctionReturn(0);
80dedccee8SHong Zhang }
81