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