1*c4762a1bSJed Brown /* This program illustrates use of parallel real FFT */ 2*c4762a1bSJed Brown static char help[]="This program illustrates the use of parallel real multi-dimensional fftw (without PETSc interface)"; 3*c4762a1bSJed Brown #include <petscmat.h> 4*c4762a1bSJed Brown #include <fftw3.h> 5*c4762a1bSJed Brown #include <fftw3-mpi.h> 6*c4762a1bSJed Brown 7*c4762a1bSJed Brown PetscInt main(PetscInt argc,char **args) 8*c4762a1bSJed Brown { 9*c4762a1bSJed Brown ptrdiff_t N0=2,N1=2,N2=2,N3=2,dim[4],N,D; 10*c4762a1bSJed Brown fftw_plan bplan,fplan; 11*c4762a1bSJed Brown fftw_complex *out; 12*c4762a1bSJed Brown double *in1,*in2; 13*c4762a1bSJed Brown ptrdiff_t alloc_local,local_n0,local_0_start; 14*c4762a1bSJed Brown ptrdiff_t local_n1,local_1_start; 15*c4762a1bSJed Brown PetscInt i,j,indx[100],n1; 16*c4762a1bSJed Brown PetscInt size,rank,n,*in,N_factor; 17*c4762a1bSJed Brown PetscScalar *data_fin,value1,one=1.0,zero=0.0; 18*c4762a1bSJed Brown PetscScalar a,*x_arr,*y_arr,*z_arr,enorm; 19*c4762a1bSJed Brown Vec fin,fout,fout1,x,y; 20*c4762a1bSJed Brown PetscRandom rnd; 21*c4762a1bSJed Brown PetscErrorCode ierr; 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 24*c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 25*c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers. Your current scalar type is complex"); 26*c4762a1bSJed Brown #endif 27*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRQ(ierr); 28*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr); 29*c4762a1bSJed Brown 30*c4762a1bSJed Brown PetscRandomCreate(PETSC_COMM_WORLD,&rnd); 31*c4762a1bSJed Brown D =4; 32*c4762a1bSJed Brown dim[0]=N0;dim[1]=N1;dim[2]=N2;dim[3]=N3/2+1; 33*c4762a1bSJed Brown 34*c4762a1bSJed Brown 35*c4762a1bSJed Brown alloc_local = fftw_mpi_local_size_transposed(D,dim,PETSC_COMM_WORLD,&local_n0,&local_0_start,&local_n1,&local_1_start); 36*c4762a1bSJed Brown 37*c4762a1bSJed Brown printf("The value alloc_local is %ld from process %d\n",alloc_local,rank); 38*c4762a1bSJed Brown printf("The value local_n0 is %ld from process %d\n",local_n0,rank); 39*c4762a1bSJed Brown printf("The value local_0_start is %ld from process %d\n",local_0_start,rank); 40*c4762a1bSJed Brown printf("The value local_n1 is %ld from process %d\n",local_n1,rank); 41*c4762a1bSJed Brown printf("The value local_1_start is %ld from process %d\n",local_1_start,rank); 42*c4762a1bSJed Brown 43*c4762a1bSJed Brown /* Allocate space for input and output arrays */ 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown in1=(double*)fftw_malloc(sizeof(double)*alloc_local*2); 46*c4762a1bSJed Brown in2=(double*)fftw_malloc(sizeof(double)*alloc_local*2); 47*c4762a1bSJed Brown out=(fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local); 48*c4762a1bSJed Brown 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown N=2*N0*N1*N2*(N3/2+1);N_factor=N0*N1*N2*N3; 51*c4762a1bSJed Brown n=2*local_n0*N1*N2*(N3/2+1);n1=local_n1*N0*2*N1*N2; 52*c4762a1bSJed Brown 53*c4762a1bSJed Brown /* printf("The value N is %d from process %d\n",N,rank); */ 54*c4762a1bSJed Brown /* printf("The value n is %d from process %d\n",n,rank); */ 55*c4762a1bSJed Brown /* printf("The value n1 is %d from process %d\n",n1,rank); */ 56*c4762a1bSJed Brown /* Creating data vector and accompanying array with VeccreateMPIWithArray */ 57*c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,N,(PetscScalar*)in1,&fin);CHKERRQ(ierr); 58*c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,N,(PetscScalar*)out,&fout);CHKERRQ(ierr); 59*c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,N,(PetscScalar*)in2,&fout1);CHKERRQ(ierr); 60*c4762a1bSJed Brown 61*c4762a1bSJed Brown /* VecGetSize(fin,&size); */ 62*c4762a1bSJed Brown /* printf("The size is %d\n",size); */ 63*c4762a1bSJed Brown 64*c4762a1bSJed Brown VecSet(fin,one); 65*c4762a1bSJed Brown /* VecAssemblyBegin(fin); */ 66*c4762a1bSJed Brown /* VecAssemblyEnd(fin); */ 67*c4762a1bSJed Brown /* VecView(fin,PETSC_VIEWER_STDOUT_WORLD); */ 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown 70*c4762a1bSJed Brown VecGetArray(fin,&x_arr); 71*c4762a1bSJed Brown VecGetArray(fout1,&z_arr); 72*c4762a1bSJed Brown VecGetArray(fout,&y_arr); 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown dim[3]=N3; 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown fplan=fftw_mpi_plan_dft_r2c(D,dim,(double*)x_arr,(fftw_complex*)y_arr,PETSC_COMM_WORLD,FFTW_ESTIMATE); 77*c4762a1bSJed Brown bplan=fftw_mpi_plan_dft_c2r(D,dim,(fftw_complex*)y_arr,(double*)z_arr,PETSC_COMM_WORLD,FFTW_ESTIMATE); 78*c4762a1bSJed Brown 79*c4762a1bSJed Brown fftw_execute(fplan); 80*c4762a1bSJed Brown fftw_execute(bplan); 81*c4762a1bSJed Brown 82*c4762a1bSJed Brown VecRestoreArray(fin,&x_arr); 83*c4762a1bSJed Brown VecRestoreArray(fout1,&z_arr); 84*c4762a1bSJed Brown VecRestoreArray(fout,&y_arr); 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown /* a = 1.0/(PetscReal)N_factor; */ 87*c4762a1bSJed Brown /* ierr = VecScale(fout1,a);CHKERRQ(ierr); */ 88*c4762a1bSJed Brown 89*c4762a1bSJed Brown VecAssemblyBegin(fout1); 90*c4762a1bSJed Brown VecAssemblyEnd(fout1); 91*c4762a1bSJed Brown 92*c4762a1bSJed Brown VecView(fout1,PETSC_VIEWER_STDOUT_WORLD); 93*c4762a1bSJed Brown 94*c4762a1bSJed Brown fftw_destroy_plan(fplan); 95*c4762a1bSJed Brown fftw_destroy_plan(bplan); 96*c4762a1bSJed Brown fftw_free(in1); ierr = VecDestroy(&fin);CHKERRQ(ierr); 97*c4762a1bSJed Brown fftw_free(out); ierr = VecDestroy(&fout);CHKERRQ(ierr); 98*c4762a1bSJed Brown fftw_free(in2); ierr = VecDestroy(&fout1);CHKERRQ(ierr); 99*c4762a1bSJed Brown 100*c4762a1bSJed Brown ierr = PetscFinalize(); 101*c4762a1bSJed Brown return ierr; 102*c4762a1bSJed Brown } 103