xref: /petsc/src/mat/impls/fft/fft.c (revision b00a91154f763f12aa55f3d53a3f2776f15f49e3)
1dedccee8SHong Zhang 
2dedccee8SHong Zhang /*
3dedccee8SHong Zhang     Provides an interface to the FFT packages.
4dedccee8SHong Zhang */
5dedccee8SHong Zhang 
6c6db04a5SJed Brown #include <../src/mat/impls/fft/fft.h>   /*I "petscmat.h" I*/
7dedccee8SHong Zhang 
8dedccee8SHong Zhang #undef __FUNCT__
9dedccee8SHong Zhang #define __FUNCT__ "MatDestroy_FFT"
10dedccee8SHong Zhang PetscErrorCode MatDestroy_FFT(Mat A)
11dedccee8SHong Zhang {
12dedccee8SHong Zhang   PetscErrorCode ierr;
13dedccee8SHong Zhang   Mat_FFT        *fft = (Mat_FFT*)A->data;
14dedccee8SHong Zhang 
15dedccee8SHong Zhang   PetscFunctionBegin;
16dedccee8SHong Zhang   if (fft->matdestroy) {
17dedccee8SHong Zhang     ierr = (fft->matdestroy)(A);CHKERRQ(ierr);
18dedccee8SHong Zhang   }
19dedccee8SHong Zhang   ierr = PetscFree(fft->dim);CHKERRQ(ierr);
20bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
21dedccee8SHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
22dedccee8SHong Zhang   PetscFunctionReturn(0);
23dedccee8SHong Zhang }
24dedccee8SHong Zhang 
25dedccee8SHong Zhang #undef __FUNCT__
26dedccee8SHong Zhang #define __FUNCT__ "MatCreateFFT"
27dedccee8SHong Zhang /*@
28dedccee8SHong Zhang       MatCreateFFT - Creates a matrix object that provides FFT via an external package
29dedccee8SHong Zhang 
30dedccee8SHong Zhang    Collective on MPI_Comm
31dedccee8SHong Zhang 
32dedccee8SHong Zhang    Input Parameter:
33dedccee8SHong Zhang +   comm - MPI communicator
34dedccee8SHong Zhang .   ndim - the ndim-dimensional transform
35dedccee8SHong Zhang .   dim - array of size ndim, dim[i] contains the vector length in the i-dimension
36dedccee8SHong Zhang -   type - package type, e.g., FFTW or FFTCU
37dedccee8SHong Zhang 
38dedccee8SHong Zhang    Output Parameter:
39dedccee8SHong Zhang .   A  - the matrix
40dedccee8SHong Zhang 
41dedccee8SHong Zhang   Options Database Keys:
42dedccee8SHong Zhang + -mat_fft_type - set FFT type
43dedccee8SHong Zhang 
44dedccee8SHong Zhang    Level: intermediate
45dedccee8SHong Zhang 
46dedccee8SHong Zhang @*/
4719fd82e9SBarry Smith PetscErrorCode MatCreateFFT(MPI_Comm comm,PetscInt ndim,const PetscInt dim[],MatType mattype,Mat *A)
48dedccee8SHong Zhang {
49dedccee8SHong Zhang   PetscErrorCode ierr;
50dedccee8SHong Zhang   PetscMPIInt    size;
51dedccee8SHong Zhang   Mat            FFT;
52dedccee8SHong Zhang   PetscInt       N,i;
53dedccee8SHong Zhang   Mat_FFT        *fft;
54dedccee8SHong Zhang 
55dedccee8SHong Zhang   PetscFunctionBegin;
56dedccee8SHong Zhang   if (ndim < 1) SETERRQ1(comm,PETSC_ERR_USER,"ndim %d must be > 0",ndim);
57dedccee8SHong Zhang   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
58dedccee8SHong Zhang 
59dedccee8SHong Zhang   ierr      = MatCreate(comm,&FFT);CHKERRQ(ierr);
60*b00a9115SJed Brown   ierr      = PetscNewLog(FFT,&fft);CHKERRQ(ierr);
61dedccee8SHong Zhang   FFT->data = (void*)fft;
62dedccee8SHong Zhang   N         = 1;
63dedccee8SHong Zhang   for (i=0; i<ndim; i++) {
64dedccee8SHong Zhang     if (dim[i] < 1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"dim[%d]=%d must be > 0",i,dim[i]);
65dedccee8SHong Zhang     N *= dim[i];
66dedccee8SHong Zhang   }
67dedccee8SHong Zhang 
68785e854fSJed Brown   ierr = PetscMalloc1(ndim,&fft->dim);CHKERRQ(ierr);
69dedccee8SHong Zhang   ierr = PetscMemcpy(fft->dim,dim,ndim*sizeof(PetscInt));CHKERRQ(ierr);
70dedccee8SHong Zhang 
71dedccee8SHong Zhang   fft->ndim = ndim;
72dedccee8SHong Zhang   fft->n    = PETSC_DECIDE;
73dedccee8SHong Zhang   fft->N    = N;
740298fd71SBarry Smith   fft->data = NULL;
75dedccee8SHong Zhang 
76dedccee8SHong Zhang   ierr = MatSetType(FFT,mattype);CHKERRQ(ierr);
7726fbe8dcSKarl Rupp 
78dedccee8SHong Zhang   FFT->ops->destroy = MatDestroy_FFT;
79dedccee8SHong Zhang 
80dedccee8SHong Zhang   /* get runtime options */
81ce94432eSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)FFT),((PetscObject)FFT)->prefix,"FFT Options","Mat");CHKERRQ(ierr);
82dedccee8SHong Zhang   PetscOptionsEnd();
83dedccee8SHong Zhang 
84dedccee8SHong Zhang   *A = FFT;
85dedccee8SHong Zhang   PetscFunctionReturn(0);
86dedccee8SHong Zhang }
87