xref: /petsc/src/mat/impls/fft/fft.c (revision d0609ced746bc51b019815ca91d747429db24893)
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 
7dedccee8SHong Zhang PetscErrorCode MatDestroy_FFT(Mat A)
8dedccee8SHong Zhang {
9dedccee8SHong Zhang   Mat_FFT *fft = (Mat_FFT*)A->data;
10dedccee8SHong Zhang 
11dedccee8SHong Zhang   PetscFunctionBegin;
129566063dSJacob Faibussowitsch   if (fft->matdestroy) PetscCall((fft->matdestroy)(A));
139566063dSJacob Faibussowitsch   PetscCall(PetscFree(fft->dim));
149566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
159566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A,NULL));
16dedccee8SHong Zhang   PetscFunctionReturn(0);
17dedccee8SHong Zhang }
18dedccee8SHong Zhang 
190a9977b2SMatthew G. Knepley /*@C
20dedccee8SHong Zhang       MatCreateFFT - Creates a matrix object that provides FFT via an external package
21dedccee8SHong Zhang 
22d083f849SBarry Smith    Collective
23dedccee8SHong Zhang 
24d8d19677SJose E. Roman    Input Parameters:
25dedccee8SHong Zhang +   comm - MPI communicator
26dedccee8SHong Zhang .   ndim - the ndim-dimensional transform
27dedccee8SHong Zhang .   dim - array of size ndim, dim[i] contains the vector length in the i-dimension
2875f45d78SBarry Smith -   type - package type, e.g., FFTW or MATSEQCUFFT
29dedccee8SHong Zhang 
30dedccee8SHong Zhang    Output Parameter:
31dedccee8SHong Zhang .   A  - the matrix
32dedccee8SHong Zhang 
33dedccee8SHong Zhang    Options Database Keys:
3475f45d78SBarry Smith .   -mat_fft_type - set FFT type fft or seqcufft
3575f45d78SBarry Smith 
3675f45d78SBarry Smith    Note: this serves as a base class for all FFT marix classes, currently MATFFTW or MATSEQCUFFT
37dedccee8SHong Zhang 
38dedccee8SHong Zhang    Level: intermediate
39dedccee8SHong Zhang 
405e4d437fSMatthew Knepley .seealso: MatCreateVecsFFTW()
41dedccee8SHong Zhang @*/
4219fd82e9SBarry Smith PetscErrorCode MatCreateFFT(MPI_Comm comm,PetscInt ndim,const PetscInt dim[],MatType mattype,Mat *A)
43dedccee8SHong Zhang {
44dedccee8SHong Zhang   PetscMPIInt    size;
45dedccee8SHong Zhang   Mat            FFT;
46dedccee8SHong Zhang   PetscInt       N,i;
47dedccee8SHong Zhang   Mat_FFT        *fft;
48dedccee8SHong Zhang 
49dedccee8SHong Zhang   PetscFunctionBegin;
50c0aa6a63SJacob Faibussowitsch   PetscValidIntPointer(dim,3);
51c0aa6a63SJacob Faibussowitsch   PetscValidPointer(A,5);
525f80ce2aSJacob Faibussowitsch   PetscCheck(ndim >= 1,comm,PETSC_ERR_USER,"ndim %" PetscInt_FMT " must be > 0",ndim);
539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
54dedccee8SHong Zhang 
559566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm,&FFT));
569566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(FFT,&fft));
57dedccee8SHong Zhang   FFT->data = (void*)fft;
58dedccee8SHong Zhang   N         = 1;
59dedccee8SHong Zhang   for (i=0; i<ndim; i++) {
605f80ce2aSJacob Faibussowitsch     PetscCheck(dim[i] >= 1,PETSC_COMM_SELF,PETSC_ERR_USER,"dim[%" PetscInt_FMT "]=%" PetscInt_FMT " must be > 0",i,dim[i]);
61dedccee8SHong Zhang     N *= dim[i];
62dedccee8SHong Zhang   }
63dedccee8SHong Zhang 
649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ndim,&fft->dim));
659566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(fft->dim,dim,ndim));
66dedccee8SHong Zhang 
67dedccee8SHong Zhang   fft->ndim = ndim;
68dedccee8SHong Zhang   fft->n    = PETSC_DECIDE;
69dedccee8SHong Zhang   fft->N    = N;
700298fd71SBarry Smith   fft->data = NULL;
71dedccee8SHong Zhang 
729566063dSJacob Faibussowitsch   PetscCall(MatSetType(FFT,mattype));
7326fbe8dcSKarl Rupp 
74dedccee8SHong Zhang   FFT->ops->destroy = MatDestroy_FFT;
75dedccee8SHong Zhang 
765f80ce2aSJacob Faibussowitsch   /* get runtime options... what options? */
77*d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)FFT);
78*d0609cedSBarry Smith   PetscOptionsEnd();
79dedccee8SHong Zhang 
80dedccee8SHong Zhang   *A = FFT;
81dedccee8SHong Zhang   PetscFunctionReturn(0);
82dedccee8SHong Zhang }
83