xref: /petsc/src/mat/impls/fft/fft.c (revision 5e4d437f8f9407a5762f30d616143e2f50e2f52e)
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   PetscErrorCode ierr;
10dedccee8SHong Zhang   Mat_FFT        *fft = (Mat_FFT*)A->data;
11dedccee8SHong Zhang 
12dedccee8SHong Zhang   PetscFunctionBegin;
13dedccee8SHong Zhang   if (fft->matdestroy) {
14dedccee8SHong Zhang     ierr = (fft->matdestroy)(A);CHKERRQ(ierr);
15dedccee8SHong Zhang   }
16dedccee8SHong Zhang   ierr = PetscFree(fft->dim);CHKERRQ(ierr);
17bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
18dedccee8SHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
19dedccee8SHong Zhang   PetscFunctionReturn(0);
20dedccee8SHong Zhang }
21dedccee8SHong Zhang 
220a9977b2SMatthew G. Knepley /*@C
23dedccee8SHong Zhang       MatCreateFFT - Creates a matrix object that provides FFT via an external package
24dedccee8SHong Zhang 
25d083f849SBarry Smith    Collective
26dedccee8SHong Zhang 
27dedccee8SHong Zhang    Input Parameter:
28dedccee8SHong Zhang +   comm - MPI communicator
29dedccee8SHong Zhang .   ndim - the ndim-dimensional transform
30dedccee8SHong Zhang .   dim - array of size ndim, dim[i] contains the vector length in the i-dimension
3175f45d78SBarry Smith -   type - package type, e.g., FFTW or MATSEQCUFFT
32dedccee8SHong Zhang 
33dedccee8SHong Zhang    Output Parameter:
34dedccee8SHong Zhang .   A  - the matrix
35dedccee8SHong Zhang 
36dedccee8SHong Zhang    Options Database Keys:
3775f45d78SBarry Smith .   -mat_fft_type - set FFT type fft or seqcufft
3875f45d78SBarry Smith 
3975f45d78SBarry Smith    Note: this serves as a base class for all FFT marix classes, currently MATFFTW or MATSEQCUFFT
40dedccee8SHong Zhang 
41dedccee8SHong Zhang    Level: intermediate
42dedccee8SHong Zhang 
43*5e4d437fSMatthew Knepley .seealso: MatCreateVecsFFTW()
44dedccee8SHong Zhang @*/
4519fd82e9SBarry Smith PetscErrorCode MatCreateFFT(MPI_Comm comm,PetscInt ndim,const PetscInt dim[],MatType mattype,Mat *A)
46dedccee8SHong Zhang {
47dedccee8SHong Zhang   PetscErrorCode ierr;
48dedccee8SHong Zhang   PetscMPIInt    size;
49dedccee8SHong Zhang   Mat            FFT;
50dedccee8SHong Zhang   PetscInt       N,i;
51dedccee8SHong Zhang   Mat_FFT        *fft;
52dedccee8SHong Zhang 
53dedccee8SHong Zhang   PetscFunctionBegin;
54dedccee8SHong Zhang   if (ndim < 1) SETERRQ1(comm,PETSC_ERR_USER,"ndim %d must be > 0",ndim);
55dedccee8SHong Zhang   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
56dedccee8SHong Zhang 
57dedccee8SHong Zhang   ierr      = MatCreate(comm,&FFT);CHKERRQ(ierr);
58b00a9115SJed Brown   ierr      = PetscNewLog(FFT,&fft);CHKERRQ(ierr);
59dedccee8SHong Zhang   FFT->data = (void*)fft;
60dedccee8SHong Zhang   N         = 1;
61dedccee8SHong Zhang   for (i=0; i<ndim; i++) {
62dedccee8SHong Zhang     if (dim[i] < 1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"dim[%d]=%d must be > 0",i,dim[i]);
63dedccee8SHong Zhang     N *= dim[i];
64dedccee8SHong Zhang   }
65dedccee8SHong Zhang 
66785e854fSJed Brown   ierr = PetscMalloc1(ndim,&fft->dim);CHKERRQ(ierr);
67580bdb30SBarry Smith   ierr = PetscArraycpy(fft->dim,dim,ndim);CHKERRQ(ierr);
68dedccee8SHong Zhang 
69dedccee8SHong Zhang   fft->ndim = ndim;
70dedccee8SHong Zhang   fft->n    = PETSC_DECIDE;
71dedccee8SHong Zhang   fft->N    = N;
720298fd71SBarry Smith   fft->data = NULL;
73dedccee8SHong Zhang 
74dedccee8SHong Zhang   ierr = MatSetType(FFT,mattype);CHKERRQ(ierr);
7526fbe8dcSKarl Rupp 
76dedccee8SHong Zhang   FFT->ops->destroy = MatDestroy_FFT;
77dedccee8SHong Zhang 
78dedccee8SHong Zhang   /* get runtime options */
79ce94432eSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)FFT),((PetscObject)FFT)->prefix,"FFT Options","Mat");CHKERRQ(ierr);
80dedccee8SHong Zhang   PetscOptionsEnd();
81dedccee8SHong Zhang 
82dedccee8SHong Zhang   *A = FFT;
83dedccee8SHong Zhang   PetscFunctionReturn(0);
84dedccee8SHong Zhang }
85