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