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 2D fft using fftw without PETSc interface"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown #include <petscmat.h> 5*c4762a1bSJed Brown #include <fftw3.h> 6*c4762a1bSJed Brown #include <fftw3-mpi.h> 7*c4762a1bSJed Brown 8*c4762a1bSJed Brown int main(int argc,char **args) 9*c4762a1bSJed Brown { 10*c4762a1bSJed Brown const ptrdiff_t N0=2056,N1=2056; 11*c4762a1bSJed Brown fftw_plan bplan,fplan; 12*c4762a1bSJed Brown fftw_complex *out; 13*c4762a1bSJed Brown double *in1,*in2; 14*c4762a1bSJed Brown ptrdiff_t alloc_local,local_n0,local_0_start; 15*c4762a1bSJed Brown ptrdiff_t local_n1,local_1_start; 16*c4762a1bSJed Brown PetscInt i,j; 17*c4762a1bSJed Brown PetscMPIInt size,rank; 18*c4762a1bSJed Brown int n,N,N_factor,NM; 19*c4762a1bSJed Brown PetscScalar one=2.0,zero=0.5; 20*c4762a1bSJed Brown PetscScalar two=4.0,three=8.0,four=16.0; 21*c4762a1bSJed Brown PetscScalar a,*x_arr,*y_arr,*z_arr; 22*c4762a1bSJed Brown PetscReal enorm; 23*c4762a1bSJed Brown Vec fin,fout,fout1; 24*c4762a1bSJed Brown Vec ini,final; 25*c4762a1bSJed Brown PetscRandom rnd; 26*c4762a1bSJed Brown PetscErrorCode ierr; 27*c4762a1bSJed Brown PetscInt *indx3,tempindx,low,*indx4,tempindx1; 28*c4762a1bSJed Brown 29*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 30*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr); 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rnd);CHKERRQ(ierr); 34*c4762a1bSJed Brown 35*c4762a1bSJed Brown alloc_local = fftw_mpi_local_size_2d_transposed(N0,N1/2+1,PETSC_COMM_WORLD,&local_n0,&local_0_start,&local_n1,&local_1_start); 36*c4762a1bSJed Brown #if defined(DEBUGGING) 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 /* printf("The value local_n0 is %ld from process %d\n",local_n0,rank); */ 43*c4762a1bSJed Brown #endif 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown /* Allocate space for input and output arrays */ 46*c4762a1bSJed Brown in1=(double*)fftw_malloc(sizeof(double)*alloc_local*2); 47*c4762a1bSJed Brown in2=(double*)fftw_malloc(sizeof(double)*alloc_local*2); 48*c4762a1bSJed Brown out=(fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local); 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown N = 2*N0*(N1/2+1); 51*c4762a1bSJed Brown N_factor = N0*N1; 52*c4762a1bSJed Brown n = 2*local_n0*(N1/2+1); 53*c4762a1bSJed Brown 54*c4762a1bSJed Brown /* printf("The value N is %d from process %d\n",N,rank); */ 55*c4762a1bSJed Brown /* printf("The value n is %d from process %d\n",n,rank); */ 56*c4762a1bSJed Brown /* printf("The value n1 is %d from process %d\n",n1,rank);*/ 57*c4762a1bSJed Brown /* Creating data vector and accompanying array with VeccreateMPIWithArray */ 58*c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,N,(PetscScalar*)in1,&fin);CHKERRQ(ierr); 59*c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,N,(PetscScalar*)out,&fout);CHKERRQ(ierr); 60*c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,N,(PetscScalar*)in2,&fout1);CHKERRQ(ierr); 61*c4762a1bSJed Brown 62*c4762a1bSJed Brown /* Set the vector with random data */ 63*c4762a1bSJed Brown ierr = VecSet(fin,zero);CHKERRQ(ierr); 64*c4762a1bSJed Brown /* for (i=0;i<N0*N1;i++) */ 65*c4762a1bSJed Brown /* { */ 66*c4762a1bSJed Brown /* VecSetValues(fin,1,&i,&one,INSERT_VALUES); */ 67*c4762a1bSJed Brown /* } */ 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown /* VecSet(fin,one); */ 70*c4762a1bSJed Brown i =0; 71*c4762a1bSJed Brown ierr = VecSetValues(fin,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 72*c4762a1bSJed Brown i =1; 73*c4762a1bSJed Brown ierr = VecSetValues(fin,1,&i,&two,INSERT_VALUES);CHKERRQ(ierr); 74*c4762a1bSJed Brown i =4; 75*c4762a1bSJed Brown ierr = VecSetValues(fin,1,&i,&three,INSERT_VALUES);CHKERRQ(ierr); 76*c4762a1bSJed Brown i =5; 77*c4762a1bSJed Brown ierr = VecSetValues(fin,1,&i,&four,INSERT_VALUES);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = VecAssemblyBegin(fin);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = VecAssemblyEnd(fin);CHKERRQ(ierr); 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown ierr = VecSet(fout,zero);CHKERRQ(ierr); 82*c4762a1bSJed Brown ierr = VecSet(fout1,zero);CHKERRQ(ierr); 83*c4762a1bSJed Brown 84*c4762a1bSJed Brown /* Get the meaningful portion of array */ 85*c4762a1bSJed Brown ierr = VecGetArray(fin,&x_arr);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = VecGetArray(fout1,&z_arr);CHKERRQ(ierr); 87*c4762a1bSJed Brown ierr = VecGetArray(fout,&y_arr);CHKERRQ(ierr); 88*c4762a1bSJed Brown 89*c4762a1bSJed Brown fplan=fftw_mpi_plan_dft_r2c_2d(N0,N1,(double*)x_arr,(fftw_complex*)y_arr,PETSC_COMM_WORLD,FFTW_ESTIMATE); 90*c4762a1bSJed Brown bplan=fftw_mpi_plan_dft_c2r_2d(N0,N1,(fftw_complex*)y_arr,(double*)z_arr,PETSC_COMM_WORLD,FFTW_ESTIMATE); 91*c4762a1bSJed Brown 92*c4762a1bSJed Brown fftw_execute(fplan); 93*c4762a1bSJed Brown fftw_execute(bplan); 94*c4762a1bSJed Brown 95*c4762a1bSJed Brown ierr = VecRestoreArray(fin,&x_arr); 96*c4762a1bSJed Brown ierr = VecRestoreArray(fout1,&z_arr); 97*c4762a1bSJed Brown ierr = VecRestoreArray(fout,&y_arr); 98*c4762a1bSJed Brown 99*c4762a1bSJed Brown /* VecView(fin,PETSC_VIEWER_STDOUT_WORLD); */ 100*c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&ini);CHKERRQ(ierr); 101*c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&final);CHKERRQ(ierr); 102*c4762a1bSJed Brown ierr = VecSetSizes(ini,local_n0*N1,N0*N1);CHKERRQ(ierr); 103*c4762a1bSJed Brown ierr = VecSetSizes(final,local_n0*N1,N0*N1);CHKERRQ(ierr); 104*c4762a1bSJed Brown ierr = VecSetFromOptions(ini);CHKERRQ(ierr); 105*c4762a1bSJed Brown ierr = VecSetFromOptions(final);CHKERRQ(ierr); 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown if (N1%2==0) { 108*c4762a1bSJed Brown NM = N1+2; 109*c4762a1bSJed Brown } else { 110*c4762a1bSJed Brown NM = N1+1; 111*c4762a1bSJed Brown } 112*c4762a1bSJed Brown /*printf("The Value of NM is %d",NM); */ 113*c4762a1bSJed Brown ierr = VecGetOwnershipRange(fin,&low,NULL); 114*c4762a1bSJed Brown /*printf("The local index is %d from %d\n",low,rank); */ 115*c4762a1bSJed Brown ierr = PetscMalloc1(local_n0*N1,&indx3); 116*c4762a1bSJed Brown ierr = PetscMalloc1(local_n0*N1,&indx4); 117*c4762a1bSJed Brown for (i=0;i<local_n0;i++) { 118*c4762a1bSJed Brown for (j=0;j<N1;j++) { 119*c4762a1bSJed Brown tempindx = i*N1 + j; 120*c4762a1bSJed Brown tempindx1 = i*NM + j; 121*c4762a1bSJed Brown 122*c4762a1bSJed Brown indx3[tempindx]=local_0_start*N1+tempindx; 123*c4762a1bSJed Brown indx4[tempindx]=low+tempindx1; 124*c4762a1bSJed Brown /* printf("index3 %d from proc %d is \n",indx3[tempindx],rank); */ 125*c4762a1bSJed Brown /* printf("index4 %d from proc %d is \n",indx4[tempindx],rank); */ 126*c4762a1bSJed Brown } 127*c4762a1bSJed Brown } 128*c4762a1bSJed Brown 129*c4762a1bSJed Brown ierr = PetscMalloc2(local_n0*N1,&x_arr,local_n0*N1,&y_arr);CHKERRQ(ierr); /* arr must be allocated for VecGetValues() */ 130*c4762a1bSJed Brown ierr = VecGetValues(fin,local_n0*N1,indx4,(PetscScalar*)x_arr);CHKERRQ(ierr); 131*c4762a1bSJed Brown ierr = VecSetValues(ini,local_n0*N1,indx3,x_arr,INSERT_VALUES);CHKERRQ(ierr); 132*c4762a1bSJed Brown 133*c4762a1bSJed Brown ierr = VecAssemblyBegin(ini);CHKERRQ(ierr); 134*c4762a1bSJed Brown ierr = VecAssemblyEnd(ini);CHKERRQ(ierr); 135*c4762a1bSJed Brown 136*c4762a1bSJed Brown ierr = VecGetValues(fout1,local_n0*N1,indx4,y_arr); 137*c4762a1bSJed Brown ierr = VecSetValues(final,local_n0*N1,indx3,y_arr,INSERT_VALUES); 138*c4762a1bSJed Brown ierr = VecAssemblyBegin(final); 139*c4762a1bSJed Brown ierr = VecAssemblyEnd(final); 140*c4762a1bSJed Brown ierr = PetscFree2(x_arr,y_arr);CHKERRQ(ierr); 141*c4762a1bSJed Brown 142*c4762a1bSJed Brown /* 143*c4762a1bSJed Brown VecScatter vecscat; 144*c4762a1bSJed Brown IS indx1,indx2; 145*c4762a1bSJed Brown for (i=0;i<N0;i++) { 146*c4762a1bSJed Brown indx = i*NM; 147*c4762a1bSJed Brown ISCreateStride(PETSC_COMM_WORLD,N1,indx,1,&indx1); 148*c4762a1bSJed Brown indx = i*N1; 149*c4762a1bSJed Brown ISCreateStride(PETSC_COMM_WORLD,N1,indx,1,&indx2); 150*c4762a1bSJed Brown VecScatterCreate(fin,indx1,ini,indx2,&vecscat); 151*c4762a1bSJed Brown VecScatterBegin(vecscat,fin,ini,INSERT_VALUES,SCATTER_FORWARD); 152*c4762a1bSJed Brown VecScatterEnd(vecscat,fin,ini,INSERT_VALUES,SCATTER_FORWARD); 153*c4762a1bSJed Brown VecScatterBegin(vecscat,fout1,final,INSERT_VALUES,SCATTER_FORWARD); 154*c4762a1bSJed Brown VecScatterEnd(vecscat,fout1,final,INSERT_VALUES,SCATTER_FORWARD); 155*c4762a1bSJed Brown } 156*c4762a1bSJed Brown */ 157*c4762a1bSJed Brown 158*c4762a1bSJed Brown a = 1.0/(PetscReal)N_factor; 159*c4762a1bSJed Brown ierr = VecScale(fout1,a);CHKERRQ(ierr); 160*c4762a1bSJed Brown ierr = VecScale(final,a);CHKERRQ(ierr); 161*c4762a1bSJed Brown 162*c4762a1bSJed Brown 163*c4762a1bSJed Brown /* VecView(ini,PETSC_VIEWER_STDOUT_WORLD); */ 164*c4762a1bSJed Brown /* VecView(final,PETSC_VIEWER_STDOUT_WORLD); */ 165*c4762a1bSJed Brown ierr = VecAXPY(final,-1.0,ini);CHKERRQ(ierr); 166*c4762a1bSJed Brown 167*c4762a1bSJed Brown ierr = VecNorm(final,NORM_1,&enorm);CHKERRQ(ierr); 168*c4762a1bSJed Brown if (enorm > 1.e-10) { 169*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," Error norm of |x - z| = %e\n",enorm);CHKERRQ(ierr); 170*c4762a1bSJed Brown } 171*c4762a1bSJed Brown 172*c4762a1bSJed Brown /* Execute fftw with function fftw_execute and destroy it after execution */ 173*c4762a1bSJed Brown fftw_destroy_plan(fplan); 174*c4762a1bSJed Brown fftw_destroy_plan(bplan); 175*c4762a1bSJed Brown fftw_free(in1); ierr = VecDestroy(&fin);CHKERRQ(ierr); 176*c4762a1bSJed Brown fftw_free(out); ierr = VecDestroy(&fout);CHKERRQ(ierr); 177*c4762a1bSJed Brown fftw_free(in2); ierr = VecDestroy(&fout1);CHKERRQ(ierr); 178*c4762a1bSJed Brown 179*c4762a1bSJed Brown ierr = VecDestroy(&ini);CHKERRQ(ierr); 180*c4762a1bSJed Brown ierr = VecDestroy(&final);CHKERRQ(ierr); 181*c4762a1bSJed Brown 182*c4762a1bSJed Brown ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr); 183*c4762a1bSJed Brown ierr = PetscFree(indx3);CHKERRQ(ierr); 184*c4762a1bSJed Brown ierr = PetscFree(indx4);CHKERRQ(ierr); 185*c4762a1bSJed Brown ierr = PetscFinalize(); 186*c4762a1bSJed Brown return ierr; 187*c4762a1bSJed Brown } 188*c4762a1bSJed Brown 189*c4762a1bSJed Brown 190*c4762a1bSJed Brown 191*c4762a1bSJed Brown 192*c4762a1bSJed Brown /*TEST 193*c4762a1bSJed Brown 194*c4762a1bSJed Brown build: 195*c4762a1bSJed Brown requires: fftw !complex 196*c4762a1bSJed Brown 197*c4762a1bSJed Brown test: 198*c4762a1bSJed Brown output_file: output/ex144.out 199*c4762a1bSJed Brown 200*c4762a1bSJed Brown test: 201*c4762a1bSJed Brown suffix: 2 202*c4762a1bSJed Brown nsize: 3 203*c4762a1bSJed Brown output_file: output/ex144.out 204*c4762a1bSJed Brown 205*c4762a1bSJed Brown TEST*/ 206