1c4762a1bSJed Brown static char help[] = "Illustrate how to use mpi FFTW and PETSc-FFTW interface \n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Usage: 5c4762a1bSJed Brown mpiexec -n <np> ./ex158 -use_FFTW_interface NO 6c4762a1bSJed Brown mpiexec -n <np> ./ex158 -use_FFTW_interface YES 7c4762a1bSJed Brown */ 8c4762a1bSJed Brown 9c4762a1bSJed Brown #include <petscmat.h> 10c4762a1bSJed Brown #include <fftw3-mpi.h> 11c4762a1bSJed Brown 12c4762a1bSJed Brown int main(int argc,char **args) 13c4762a1bSJed Brown { 14c4762a1bSJed Brown PetscErrorCode ierr; 15c4762a1bSJed Brown PetscMPIInt rank,size; 16c4762a1bSJed Brown PetscInt N0=50,N1=20,N=N0*N1; 17c4762a1bSJed Brown PetscRandom rdm; 18c4762a1bSJed Brown PetscScalar a; 19c4762a1bSJed Brown PetscReal enorm; 20c4762a1bSJed Brown Vec x,y,z; 21c4762a1bSJed Brown PetscBool view=PETSC_FALSE,use_interface=PETSC_TRUE; 22c4762a1bSJed Brown 23c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 24c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 25c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers. Your current scalar type is complex"); 26c4762a1bSJed Brown #endif 27c4762a1bSJed Brown 28c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex158");CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = PetscOptionsBool("-use_FFTW_interface", "Use PETSc-FFTW interface", "ex158",use_interface, &use_interface, NULL);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 31c4762a1bSJed Brown 32ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRMPI(ierr); 33ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRMPI(ierr); 34c4762a1bSJed Brown 35c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD, &rdm);CHKERRQ(ierr); 36c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr); 37c4762a1bSJed Brown 38c4762a1bSJed Brown if (!use_interface) { 39c4762a1bSJed Brown /* Use mpi FFTW without PETSc-FFTW interface, 2D case only */ 40c4762a1bSJed Brown /*---------------------------------------------------------*/ 41c4762a1bSJed Brown fftw_plan fplan,bplan; 42c4762a1bSJed Brown fftw_complex *data_in,*data_out,*data_out2; 43c4762a1bSJed Brown ptrdiff_t alloc_local,local_n0,local_0_start; 44c4762a1bSJed Brown 45dd400576SPatrick Sanan if (rank == 0) printf("Use FFTW without PETSc-FFTW interface\n"); 46c4762a1bSJed Brown fftw_mpi_init(); 47c4762a1bSJed Brown N = N0*N1; 48c4762a1bSJed Brown alloc_local = fftw_mpi_local_size_2d(N0,N1,PETSC_COMM_WORLD,&local_n0,&local_0_start); 49c4762a1bSJed Brown 50c4762a1bSJed Brown data_in = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local); 51c4762a1bSJed Brown data_out = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local); 52c4762a1bSJed Brown data_out2 = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local); 53c4762a1bSJed Brown 54c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_in,&x);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) x, "Real Space vector");CHKERRQ(ierr); 56c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out,&y);CHKERRQ(ierr); 57c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) y, "Frequency space vector");CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out2,&z);CHKERRQ(ierr); 59c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) z, "Reconstructed vector");CHKERRQ(ierr); 60c4762a1bSJed Brown 61c4762a1bSJed Brown fplan = fftw_mpi_plan_dft_2d(N0,N1,data_in,data_out,PETSC_COMM_WORLD,FFTW_FORWARD,FFTW_ESTIMATE); 62c4762a1bSJed Brown bplan = fftw_mpi_plan_dft_2d(N0,N1,data_out,data_out2,PETSC_COMM_WORLD,FFTW_BACKWARD,FFTW_ESTIMATE); 63c4762a1bSJed Brown 64c4762a1bSJed Brown ierr = VecSetRandom(x, rdm);CHKERRQ(ierr); 65c4762a1bSJed Brown if (view) {ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 66c4762a1bSJed Brown 67c4762a1bSJed Brown fftw_execute(fplan); 68c4762a1bSJed Brown if (view) {ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 69c4762a1bSJed Brown 70c4762a1bSJed Brown fftw_execute(bplan); 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */ 73c4762a1bSJed Brown a = 1.0/(PetscReal)N; 74c4762a1bSJed Brown ierr = VecScale(z,a);CHKERRQ(ierr); 75c4762a1bSJed Brown if (view) {ierr = VecView(z, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 76c4762a1bSJed Brown ierr = VecAXPY(z,-1.0,x);CHKERRQ(ierr); 77c4762a1bSJed Brown ierr = VecNorm(z,NORM_1,&enorm);CHKERRQ(ierr); 78c4762a1bSJed Brown if (enorm > 1.e-11) { 79c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," Error norm of |x - z| %g\n",(double)enorm);CHKERRQ(ierr); 80c4762a1bSJed Brown } 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* Free spaces */ 83c4762a1bSJed Brown fftw_destroy_plan(fplan); 84c4762a1bSJed Brown fftw_destroy_plan(bplan); 85c4762a1bSJed Brown fftw_free(data_in); ierr = VecDestroy(&x);CHKERRQ(ierr); 86c4762a1bSJed Brown fftw_free(data_out); ierr = VecDestroy(&y);CHKERRQ(ierr); 87c4762a1bSJed Brown fftw_free(data_out2);ierr = VecDestroy(&z);CHKERRQ(ierr); 88c4762a1bSJed Brown 89c4762a1bSJed Brown } else { 90c4762a1bSJed Brown /* Use PETSc-FFTW interface */ 91c4762a1bSJed Brown /*-------------------------------------------*/ 92c4762a1bSJed Brown PetscInt i,*dim,k,DIM; 93c4762a1bSJed Brown Mat A; 94c4762a1bSJed Brown Vec input,output; 95c4762a1bSJed Brown 96c4762a1bSJed Brown N=30; 97c4762a1bSJed Brown for (i=2; i<3; i++) { /* (i=3,4: -- error in VecScatterPetscToFFTW(A,input,x); */ 98c4762a1bSJed Brown DIM = i; 99c4762a1bSJed Brown ierr = PetscMalloc1(i,&dim);CHKERRQ(ierr); 100c4762a1bSJed Brown for (k=0; k<i; k++) { 101c4762a1bSJed Brown dim[k]=30; 102c4762a1bSJed Brown } 103c4762a1bSJed Brown N *= dim[i-1]; 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* Create FFTW object */ 106dd400576SPatrick Sanan if (rank == 0) { 107c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Use PETSc-FFTW interface...%d-DIM:%d \n",DIM,N);CHKERRQ(ierr); 108c4762a1bSJed Brown } 109c4762a1bSJed Brown ierr = MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A);CHKERRQ(ierr); 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* Create FFTW vectors that are compatible with parallel layout of A */ 112c4762a1bSJed Brown ierr = MatCreateVecsFFTW(A,&x,&y,&z);CHKERRQ(ierr); 113c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) x, "Real space vector");CHKERRQ(ierr); 114c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) y, "Frequency space vector");CHKERRQ(ierr); 115c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) z, "Reconstructed vector");CHKERRQ(ierr); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* Create and set PETSc vector */ 118c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&input);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = VecSetSizes(input,PETSC_DECIDE,N);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = VecSetFromOptions(input);CHKERRQ(ierr); 121c4762a1bSJed Brown ierr = VecSetRandom(input,rdm);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = VecDuplicate(input,&output);CHKERRQ(ierr); 123c4762a1bSJed Brown if (view) {ierr = VecView(input,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* Vector input is copied to another vector x using VecScatterPetscToFFTW. This is because the user data 126c4762a1bSJed Brown can have any parallel layout. But FFTW requires special parallel layout of the data. Hence the original 127c4762a1bSJed Brown data which is in the vector "input" here, needs to be copied to a vector x, which has the correct parallel 128c4762a1bSJed Brown layout for FFTW. Also, during parallel real transform, this pads extra zeros automatically 129c4762a1bSJed Brown at the end of last dimension. This padding is required by FFTW to perform parallel real D.F.T. */ 130c4762a1bSJed Brown ierr = VecScatterPetscToFFTW(A,input,x);CHKERRQ(ierr);/* buggy for dim = 3, 4... */ 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* Apply FFTW_FORWARD and FFTW_BACKWARD */ 133c4762a1bSJed Brown ierr = MatMult(A,x,y);CHKERRQ(ierr); 134c4762a1bSJed Brown if (view) {ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 135c4762a1bSJed Brown ierr = MatMultTranspose(A,y,z);CHKERRQ(ierr); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* Output from Backward DFT needs to be modified to obtain user readable data the routine VecScatterFFTWToPetsc 138c4762a1bSJed Brown performs the job. In some sense this is the reverse operation of VecScatterPetscToFFTW. This routine gets rid of 139c4762a1bSJed Brown the extra spaces that were artificially padded to perform real parallel transform. */ 140c4762a1bSJed Brown ierr = VecScatterFFTWToPetsc(A,z,output);CHKERRQ(ierr); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */ 143c4762a1bSJed Brown a = 1.0/(PetscReal)N; 144c4762a1bSJed Brown ierr = VecScale(output,a);CHKERRQ(ierr); 145c4762a1bSJed Brown if (view) {ierr = VecView(output,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 146c4762a1bSJed Brown ierr = VecAXPY(output,-1.0,input);CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = VecNorm(output,NORM_1,&enorm);CHKERRQ(ierr); 148dd400576SPatrick Sanan if (enorm > 1.e-09 && rank == 0) { 149c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," Error norm of |x - z| %e\n",enorm);CHKERRQ(ierr); 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* Free spaces */ 153c4762a1bSJed Brown ierr = PetscFree(dim);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = VecDestroy(&input);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = VecDestroy(&output);CHKERRQ(ierr); 156c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = VecDestroy(&y);CHKERRQ(ierr); 158c4762a1bSJed Brown ierr = VecDestroy(&z);CHKERRQ(ierr); 159c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 160c4762a1bSJed Brown } 161c4762a1bSJed Brown } 162c4762a1bSJed Brown ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = PetscFinalize(); 164c4762a1bSJed Brown return ierr; 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167c4762a1bSJed Brown /*TEST 168c4762a1bSJed Brown 169c4762a1bSJed Brown build: 170*0cf2e031SBarry Smith requires: !mpiuni fftw !complex 171c4762a1bSJed Brown 172c4762a1bSJed Brown test: 173c4762a1bSJed Brown output_file: output/ex158.out 174c4762a1bSJed Brown 175c4762a1bSJed Brown test: 176c4762a1bSJed Brown suffix: 2 177c4762a1bSJed Brown nsize: 3 178c4762a1bSJed Brown 179c4762a1bSJed Brown TEST*/ 180