1*c4762a1bSJed Brown /* This program illustrates use of paralllel real FFT*/ 2*c4762a1bSJed Brown static char help[]="This program illustrates the use of parallel real 3D 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=256,N1=256,N2=256,N3=2,dim[4]; 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,n1; 16*c4762a1bSJed Brown PetscInt size,rank,n,N,*in,N_factor,NM; 17*c4762a1bSJed Brown PetscScalar *data_fin,value1,one=1.57,zero=0.0; 18*c4762a1bSJed Brown PetscScalar a,*x_arr,*y_arr,*z_arr,enorm; 19*c4762a1bSJed Brown Vec fin,fout,fout1,ini,final; 20*c4762a1bSJed Brown PetscRandom rnd; 21*c4762a1bSJed Brown PetscErrorCode ierr; 22*c4762a1bSJed Brown VecScatter vecscat,vecscat1; 23*c4762a1bSJed Brown IS indx1,indx2; 24*c4762a1bSJed Brown PetscInt *indx3,k,l,*indx4; 25*c4762a1bSJed Brown PetscInt low,tempindx,tempindx1; 26*c4762a1bSJed Brown 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 29*c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 30*c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers. Your current scalar type is complex"); 31*c4762a1bSJed Brown #endif 32*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRQ(ierr); 33*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr); 34*c4762a1bSJed Brown 35*c4762a1bSJed Brown PetscRandomCreate(PETSC_COMM_WORLD,&rnd); 36*c4762a1bSJed Brown 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown alloc_local = fftw_mpi_local_size_3d_transposed(N0,N1,N2/2+1,PETSC_COMM_WORLD,&local_n0,&local_0_start,&local_n1,&local_1_start); 39*c4762a1bSJed Brown 40*c4762a1bSJed Brown /* printf("The value alloc_local is %ld from process %d\n",alloc_local,rank); */ 41*c4762a1bSJed Brown printf("The value local_n0 is %ld from process %d\n",local_n0,rank); 42*c4762a1bSJed Brown /* printf("The value local_0_start is %ld from process %d\n",local_0_start,rank);*/ 43*c4762a1bSJed Brown /* printf("The value local_n1 is %ld from process %d\n",local_n1,rank); */ 44*c4762a1bSJed Brown /* printf("The value local_1_start is %ld from process %d\n",local_1_start,rank);*/ 45*c4762a1bSJed Brown 46*c4762a1bSJed Brown /* Allocate space for input and output arrays */ 47*c4762a1bSJed Brown 48*c4762a1bSJed Brown in1=(double*)fftw_malloc(sizeof(double)*alloc_local*2); 49*c4762a1bSJed Brown in2=(double*)fftw_malloc(sizeof(double)*alloc_local*2); 50*c4762a1bSJed Brown out=(fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local); 51*c4762a1bSJed Brown 52*c4762a1bSJed Brown 53*c4762a1bSJed Brown N=2*N0*N1*(N2/2+1);N_factor=N0*N1*N2; 54*c4762a1bSJed Brown n=2*local_n0*N1*(N2/2+1);n1=local_n1*N0*2*N1; 55*c4762a1bSJed Brown 56*c4762a1bSJed Brown /* printf("The value N is %d from process %d\n",N,rank); */ 57*c4762a1bSJed Brown /* printf("The value n is %d from process %d\n",n,rank); */ 58*c4762a1bSJed Brown /* printf("The value n1 is %d from process %d\n",n1,rank); */ 59*c4762a1bSJed Brown /* Creating data vector and accompanying array with VeccreateMPIWithArray */ 60*c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,N,(PetscScalar*)in1,&fin);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,N,(PetscScalar*)out,&fout);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,N,(PetscScalar*)in2,&fout1);CHKERRQ(ierr); 63*c4762a1bSJed Brown 64*c4762a1bSJed Brown /* VecGetSize(fin,&size); */ 65*c4762a1bSJed Brown /* printf("The size is %d\n",size); */ 66*c4762a1bSJed Brown 67*c4762a1bSJed Brown VecSet(fin,one); 68*c4762a1bSJed Brown VecSet(fout,zero); 69*c4762a1bSJed Brown VecSet(fout1,zero); 70*c4762a1bSJed Brown 71*c4762a1bSJed Brown VecAssemblyBegin(fin); 72*c4762a1bSJed Brown VecAssemblyEnd(fin); 73*c4762a1bSJed Brown /* VecView(fin,PETSC_VIEWER_STDOUT_WORLD); */ 74*c4762a1bSJed Brown 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown VecGetArray(fin,&x_arr); 77*c4762a1bSJed Brown VecGetArray(fout1,&z_arr); 78*c4762a1bSJed Brown VecGetArray(fout,&y_arr); 79*c4762a1bSJed Brown 80*c4762a1bSJed Brown fplan=fftw_mpi_plan_dft_r2c_3d(N0,N1,N2,(double*)x_arr,(fftw_complex*)y_arr,PETSC_COMM_WORLD,FFTW_ESTIMATE); 81*c4762a1bSJed Brown bplan=fftw_mpi_plan_dft_c2r_3d(N0,N1,N2,(fftw_complex*)y_arr,(double*)z_arr,PETSC_COMM_WORLD,FFTW_ESTIMATE); 82*c4762a1bSJed Brown 83*c4762a1bSJed Brown fftw_execute(fplan); 84*c4762a1bSJed Brown fftw_execute(bplan); 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown VecRestoreArray(fin,&x_arr); 87*c4762a1bSJed Brown VecRestoreArray(fout1,&z_arr); 88*c4762a1bSJed Brown VecRestoreArray(fout,&y_arr); 89*c4762a1bSJed Brown 90*c4762a1bSJed Brown 91*c4762a1bSJed Brown /* a = 1.0/(PetscReal)N_factor; */ 92*c4762a1bSJed Brown /* ierr = VecScale(fout1,a);CHKERRQ(ierr); */ 93*c4762a1bSJed Brown VecCreate(PETSC_COMM_WORLD,&ini); 94*c4762a1bSJed Brown VecCreate(PETSC_COMM_WORLD,&final); 95*c4762a1bSJed Brown VecSetSizes(ini,local_n0*N1*N2,N_factor); 96*c4762a1bSJed Brown VecSetSizes(final,local_n0*N1*N2,N_factor); 97*c4762a1bSJed Brown /* VecSetSizes(ini,PETSC_DECIDE,N_factor); */ 98*c4762a1bSJed Brown /* VecSetSizes(final,PETSC_DECIDE,N_factor); */ 99*c4762a1bSJed Brown VecSetFromOptions(ini); 100*c4762a1bSJed Brown VecSetFromOptions(final); 101*c4762a1bSJed Brown 102*c4762a1bSJed Brown if (N2%2==0) NM=N2+2; 103*c4762a1bSJed Brown else NM=N2+1; 104*c4762a1bSJed Brown 105*c4762a1bSJed Brown ierr = VecGetOwnershipRange(fin,&low,NULL); 106*c4762a1bSJed Brown printf("The local index is %d from %d\n",low,rank); 107*c4762a1bSJed Brown ierr = PetscMalloc1(local_n0*N1*N2,&indx3); 108*c4762a1bSJed Brown ierr = PetscMalloc1(local_n0*N1*N2,&indx4); 109*c4762a1bSJed Brown for (i=0; i<local_n0; i++) { 110*c4762a1bSJed Brown for (j=0;j<N1;j++) { 111*c4762a1bSJed Brown for (k=0;k<N2;k++) { 112*c4762a1bSJed Brown tempindx = i*N1*N2 + j*N2 + k; 113*c4762a1bSJed Brown tempindx1 = i*N1*NM + j*NM + k; 114*c4762a1bSJed Brown 115*c4762a1bSJed Brown indx3[tempindx]=local_0_start*N1*N2+tempindx; 116*c4762a1bSJed Brown indx4[tempindx]=low+tempindx1; 117*c4762a1bSJed Brown } 118*c4762a1bSJed Brown /* printf("index3 %d from proc %d is \n",indx3[tempindx],rank); */ 119*c4762a1bSJed Brown /* printf("index4 %d from proc %d is \n",indx4[tempindx],rank); */ 120*c4762a1bSJed Brown } 121*c4762a1bSJed Brown } 122*c4762a1bSJed Brown VecGetValues(fin,local_n0*N1*N2,indx4,x_arr); 123*c4762a1bSJed Brown VecSetValues(ini,local_n0*N1*N2,indx3,x_arr,INSERT_VALUES); 124*c4762a1bSJed Brown VecAssemblyBegin(ini); 125*c4762a1bSJed Brown VecAssemblyEnd(ini); 126*c4762a1bSJed Brown 127*c4762a1bSJed Brown VecGetValues(fout1,local_n0*N1*N2,indx4,y_arr); 128*c4762a1bSJed Brown VecSetValues(final,local_n0*N1*N2,indx3,y_arr,INSERT_VALUES); 129*c4762a1bSJed Brown VecAssemblyBegin(final); 130*c4762a1bSJed Brown VecAssemblyEnd(final); 131*c4762a1bSJed Brown 132*c4762a1bSJed Brown printf("The local index value is %ld from %d",local_n0*N1*N2,rank); 133*c4762a1bSJed Brown /* 134*c4762a1bSJed Brown for (i=0;i<N0;i++) { 135*c4762a1bSJed Brown for (j=0;j<N1;j++) { 136*c4762a1bSJed Brown indx=i*N1*NM+j*NM; 137*c4762a1bSJed Brown ISCreateStride(PETSC_COMM_WORLD,N2,indx,1,&indx1); 138*c4762a1bSJed Brown indx=i*N1*N2+j*N2; 139*c4762a1bSJed Brown ISCreateStride(PETSC_COMM_WORLD,N2,indx,1,&indx2); 140*c4762a1bSJed Brown VecScatterCreate(fin,indx1,ini,indx2,&vecscat); 141*c4762a1bSJed Brown VecScatterBegin(vecscat,fin,ini,INSERT_VALUES,SCATTER_FORWARD); 142*c4762a1bSJed Brown VecScatterEnd(vecscat,fin,ini,INSERT_VALUES,SCATTER_FORWARD); 143*c4762a1bSJed Brown VecScatterCreate(fout1,indx1,final,indx2,&vecscat1); 144*c4762a1bSJed Brown VecScatterBegin(vecscat1,fout1,final,INSERT_VALUES,SCATTER_FORWARD); 145*c4762a1bSJed Brown VecScatterEnd(vecscat1,fout1,final,INSERT_VALUES,SCATTER_FORWARD); 146*c4762a1bSJed Brown } 147*c4762a1bSJed Brown } 148*c4762a1bSJed Brown */ 149*c4762a1bSJed Brown a = 1.0/(PetscReal)N_factor; 150*c4762a1bSJed Brown ierr = VecScale(fout1,a);CHKERRQ(ierr); 151*c4762a1bSJed Brown ierr = VecScale(final,a);CHKERRQ(ierr); 152*c4762a1bSJed Brown 153*c4762a1bSJed Brown VecAssemblyBegin(ini); 154*c4762a1bSJed Brown VecAssemblyEnd(ini); 155*c4762a1bSJed Brown 156*c4762a1bSJed Brown VecAssemblyBegin(final); 157*c4762a1bSJed Brown VecAssemblyEnd(final); 158*c4762a1bSJed Brown 159*c4762a1bSJed Brown /* VecView(final,PETSC_VIEWER_STDOUT_WORLD); */ 160*c4762a1bSJed Brown ierr = VecAXPY(final,-1.0,ini);CHKERRQ(ierr); 161*c4762a1bSJed Brown ierr = VecNorm(final,NORM_1,&enorm);CHKERRQ(ierr); 162*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," Error norm of |x - z| = %e\n",enorm);CHKERRQ(ierr); 163*c4762a1bSJed Brown fftw_destroy_plan(fplan); 164*c4762a1bSJed Brown fftw_destroy_plan(bplan); 165*c4762a1bSJed Brown fftw_free(in1); ierr = VecDestroy(&fin);CHKERRQ(ierr); 166*c4762a1bSJed Brown fftw_free(out); ierr = VecDestroy(&fout);CHKERRQ(ierr); 167*c4762a1bSJed Brown fftw_free(in2); ierr = VecDestroy(&fout1);CHKERRQ(ierr); 168*c4762a1bSJed Brown 169*c4762a1bSJed Brown ierr = PetscFinalize(); 170*c4762a1bSJed Brown return ierr; 171*c4762a1bSJed Brown } 172