1c4762a1bSJed Brown static char help[] = "Illustrate how to use mpi FFTW and PETSc-FFTW interface \n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Compiling the code: 5c4762a1bSJed Brown This code uses the complex numbers version of PETSc, so configure 6c4762a1bSJed Brown must be run to enable this 7c4762a1bSJed Brown 8c4762a1bSJed Brown Usage: 9c4762a1bSJed Brown mpiexec -n <np> ./ex143 -use_FFTW_interface NO 10c4762a1bSJed Brown mpiexec -n <np> ./ex143 -use_FFTW_interface YES 11c4762a1bSJed Brown */ 12c4762a1bSJed Brown 13c4762a1bSJed Brown #include <petscmat.h> 14c4762a1bSJed Brown #include <fftw3-mpi.h> 15c4762a1bSJed Brown 16c4762a1bSJed Brown int main(int argc,char **args) 17c4762a1bSJed Brown { 18c4762a1bSJed Brown PetscErrorCode ierr; 19c4762a1bSJed Brown PetscMPIInt rank,size; 20c4762a1bSJed Brown PetscInt N0=50,N1=20,N=N0*N1,DIM; 21c4762a1bSJed Brown PetscRandom rdm; 22c4762a1bSJed Brown PetscScalar a; 23c4762a1bSJed Brown PetscReal enorm; 24c4762a1bSJed Brown Vec x,y,z; 25c4762a1bSJed Brown PetscBool view=PETSC_FALSE,use_interface=PETSC_TRUE; 26c4762a1bSJed Brown 27c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 28c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex143");CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = PetscOptionsBool("-vec_view draw", "View the vectors", "ex143", view, &view, NULL);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = PetscOptionsBool("-use_FFTW_interface", "Use PETSc-FFTW interface", "ex143",use_interface, &use_interface, NULL);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 32c4762a1bSJed Brown 33c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-use_FFTW_interface",&use_interface,NULL);CHKERRQ(ierr); 34ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRMPI(ierr); 35ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRMPI(ierr); 36c4762a1bSJed Brown 37c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD, &rdm);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr); 39c4762a1bSJed Brown 40c4762a1bSJed Brown if (!use_interface) { 41c4762a1bSJed Brown /* Use mpi FFTW without PETSc-FFTW interface, 2D case only */ 42c4762a1bSJed Brown /*---------------------------------------------------------*/ 43c4762a1bSJed Brown fftw_plan fplan,bplan; 44c4762a1bSJed Brown fftw_complex *data_in,*data_out,*data_out2; 45c4762a1bSJed Brown ptrdiff_t alloc_local,local_n0,local_0_start; 46c4762a1bSJed Brown 47c4762a1bSJed Brown DIM = 2; 48dd400576SPatrick Sanan if (rank == 0) { 49c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Use FFTW without PETSc-FFTW interface, DIM %D\n",DIM);CHKERRQ(ierr); 50c4762a1bSJed Brown } 51c4762a1bSJed Brown fftw_mpi_init(); 52c4762a1bSJed Brown N = N0*N1; 53c4762a1bSJed Brown alloc_local = fftw_mpi_local_size_2d(N0,N1,PETSC_COMM_WORLD,&local_n0,&local_0_start); 54c4762a1bSJed Brown 55c4762a1bSJed Brown data_in = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local); 56c4762a1bSJed Brown data_out = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local); 57c4762a1bSJed Brown data_out2 = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local); 58c4762a1bSJed Brown 59c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_in,&x);CHKERRQ(ierr); 60c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) x, "Real Space vector");CHKERRQ(ierr); 61c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out,&y);CHKERRQ(ierr); 62c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) y, "Frequency space vector");CHKERRQ(ierr); 63c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out2,&z);CHKERRQ(ierr); 64c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) z, "Reconstructed vector");CHKERRQ(ierr); 65c4762a1bSJed Brown 66c4762a1bSJed Brown fplan = fftw_mpi_plan_dft_2d(N0,N1,data_in,data_out,PETSC_COMM_WORLD,FFTW_FORWARD,FFTW_ESTIMATE); 67c4762a1bSJed Brown bplan = fftw_mpi_plan_dft_2d(N0,N1,data_out,data_out2,PETSC_COMM_WORLD,FFTW_BACKWARD,FFTW_ESTIMATE); 68c4762a1bSJed Brown 69c4762a1bSJed Brown ierr = VecSetRandom(x, rdm);CHKERRQ(ierr); 70c4762a1bSJed Brown if (view) {ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 71c4762a1bSJed Brown 72c4762a1bSJed Brown fftw_execute(fplan); 73c4762a1bSJed Brown if (view) {ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 74c4762a1bSJed Brown 75c4762a1bSJed Brown fftw_execute(bplan); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */ 78c4762a1bSJed Brown a = 1.0/(PetscReal)N; 79c4762a1bSJed Brown ierr = VecScale(z,a);CHKERRQ(ierr); 80c4762a1bSJed Brown if (view) {ierr = VecView(z, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 81c4762a1bSJed Brown ierr = VecAXPY(z,-1.0,x);CHKERRQ(ierr); 82c4762a1bSJed Brown ierr = VecNorm(z,NORM_1,&enorm);CHKERRQ(ierr); 83dd400576SPatrick Sanan if (enorm > 1.e-11 && rank == 0) { 84c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," Error norm of |x - z| %g\n",(double)enorm);CHKERRQ(ierr); 85c4762a1bSJed Brown } 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* Free spaces */ 88c4762a1bSJed Brown fftw_destroy_plan(fplan); 89c4762a1bSJed Brown fftw_destroy_plan(bplan); 90c4762a1bSJed Brown fftw_free(data_in); ierr = VecDestroy(&x);CHKERRQ(ierr); 91c4762a1bSJed Brown fftw_free(data_out); ierr = VecDestroy(&y);CHKERRQ(ierr); 92c4762a1bSJed Brown fftw_free(data_out2);ierr = VecDestroy(&z);CHKERRQ(ierr); 93c4762a1bSJed Brown 94c4762a1bSJed Brown } else { 95c4762a1bSJed Brown /* Use PETSc-FFTW interface */ 96c4762a1bSJed Brown /*-------------------------------------------*/ 97c4762a1bSJed Brown PetscInt i,*dim,k; 98c4762a1bSJed Brown Mat A; 99c4762a1bSJed Brown 100c4762a1bSJed Brown N=1; 101c4762a1bSJed Brown for (i=1; i<5; i++) { 102c4762a1bSJed Brown DIM = i; 103c4762a1bSJed Brown ierr = PetscMalloc1(i,&dim);CHKERRQ(ierr); 104c4762a1bSJed Brown for (k=0; k<i; k++) { 105c4762a1bSJed Brown dim[k]=30; 106c4762a1bSJed Brown } 107c4762a1bSJed Brown N *= dim[i-1]; 108c4762a1bSJed Brown 109c4762a1bSJed Brown /* Create FFTW object */ 110dd400576SPatrick Sanan if (rank == 0) printf("Use PETSc-FFTW interface...%d-DIM: %d\n",(int)DIM,(int)N); 111c4762a1bSJed Brown 112c4762a1bSJed Brown ierr = MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A);CHKERRQ(ierr); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* Create vectors that are compatible with parallel layout of A - must call MatCreateVecs()! */ 115c4762a1bSJed Brown 116c4762a1bSJed Brown ierr = MatCreateVecsFFTW(A,&x,&y,&z);CHKERRQ(ierr); 117c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) x, "Real space vector");CHKERRQ(ierr); 118c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) y, "Frequency space vector");CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) z, "Reconstructed vector");CHKERRQ(ierr); 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* Set values of space vector x */ 122c4762a1bSJed Brown ierr = VecSetRandom(x,rdm);CHKERRQ(ierr); 123c4762a1bSJed Brown 124c4762a1bSJed Brown if (view) {ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* Apply FFTW_FORWARD and FFTW_BACKWARD */ 127c4762a1bSJed Brown ierr = MatMult(A,x,y);CHKERRQ(ierr); 128c4762a1bSJed Brown if (view) {ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 129c4762a1bSJed Brown 130c4762a1bSJed Brown ierr = MatMultTranspose(A,y,z);CHKERRQ(ierr); 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */ 133c4762a1bSJed Brown a = 1.0/(PetscReal)N; 134c4762a1bSJed Brown ierr = VecScale(z,a);CHKERRQ(ierr); 135c4762a1bSJed Brown if (view) {ierr = VecView(z,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 136c4762a1bSJed Brown ierr = VecAXPY(z,-1.0,x);CHKERRQ(ierr); 137c4762a1bSJed Brown ierr = VecNorm(z,NORM_1,&enorm);CHKERRQ(ierr); 138dd400576SPatrick Sanan if (enorm > 1.e-9 && rank == 0) { 139c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," Error norm of |x - z| %e\n",enorm);CHKERRQ(ierr); 140c4762a1bSJed Brown } 141c4762a1bSJed Brown 142c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = VecDestroy(&y);CHKERRQ(ierr); 144c4762a1bSJed Brown ierr = VecDestroy(&z);CHKERRQ(ierr); 145c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 146c4762a1bSJed Brown 147c4762a1bSJed Brown ierr = PetscFree(dim);CHKERRQ(ierr); 148c4762a1bSJed Brown } 149c4762a1bSJed Brown } 150c4762a1bSJed Brown 151c4762a1bSJed Brown ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr); 152c4762a1bSJed Brown ierr = PetscFinalize(); 153c4762a1bSJed Brown return ierr; 154c4762a1bSJed Brown } 155c4762a1bSJed Brown 156c4762a1bSJed Brown /*TEST 157c4762a1bSJed Brown 158c4762a1bSJed Brown build: 159*0cf2e031SBarry Smith requires: !mpiuni fftw complex 160c4762a1bSJed Brown 161c4762a1bSJed Brown test: 162c4762a1bSJed Brown output_file: output/ex143.out 163c4762a1bSJed Brown 164c4762a1bSJed Brown test: 165c4762a1bSJed Brown suffix: 2 166c4762a1bSJed Brown nsize: 3 167c4762a1bSJed Brown output_file: output/ex143.out 168c4762a1bSJed Brown 169c4762a1bSJed Brown TEST*/ 170