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 PetscMPIInt rank,size; 19c4762a1bSJed Brown PetscInt N0=50,N1=20,N=N0*N1,DIM; 20c4762a1bSJed Brown PetscRandom rdm; 21c4762a1bSJed Brown PetscScalar a; 22c4762a1bSJed Brown PetscReal enorm; 23c4762a1bSJed Brown Vec x,y,z; 24c4762a1bSJed Brown PetscBool view=PETSC_FALSE,use_interface=PETSC_TRUE; 25c4762a1bSJed Brown 269566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 27*d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex143"); 289566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-vec_view draw", "View the vectors", "ex143", view, &view, NULL)); 299566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_FFTW_interface", "Use PETSc-FFTW interface", "ex143",use_interface, &use_interface, NULL)); 30*d0609cedSBarry Smith PetscOptionsEnd(); 31c4762a1bSJed Brown 329566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-use_FFTW_interface",&use_interface,NULL)); 339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 35c4762a1bSJed Brown 369566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm)); 379566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 38c4762a1bSJed Brown 39c4762a1bSJed Brown if (!use_interface) { 40c4762a1bSJed Brown /* Use mpi FFTW without PETSc-FFTW interface, 2D case only */ 41c4762a1bSJed Brown /*---------------------------------------------------------*/ 42c4762a1bSJed Brown fftw_plan fplan,bplan; 43c4762a1bSJed Brown fftw_complex *data_in,*data_out,*data_out2; 44c4762a1bSJed Brown ptrdiff_t alloc_local,local_n0,local_0_start; 45c4762a1bSJed Brown 46c4762a1bSJed Brown DIM = 2; 47dd400576SPatrick Sanan if (rank == 0) { 489566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Use FFTW without PETSc-FFTW interface, DIM %" PetscInt_FMT "\n",DIM)); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown fftw_mpi_init(); 51c4762a1bSJed Brown N = N0*N1; 52c4762a1bSJed Brown alloc_local = fftw_mpi_local_size_2d(N0,N1,PETSC_COMM_WORLD,&local_n0,&local_0_start); 53c4762a1bSJed Brown 54c4762a1bSJed Brown data_in = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local); 55c4762a1bSJed Brown data_out = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local); 56c4762a1bSJed Brown data_out2 = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local); 57c4762a1bSJed Brown 589566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_in,&x)); 599566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) x, "Real Space vector")); 609566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out,&y)); 619566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) y, "Frequency space vector")); 629566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out2,&z)); 639566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) z, "Reconstructed vector")); 64c4762a1bSJed Brown 65c4762a1bSJed Brown fplan = fftw_mpi_plan_dft_2d(N0,N1,data_in,data_out,PETSC_COMM_WORLD,FFTW_FORWARD,FFTW_ESTIMATE); 66c4762a1bSJed Brown bplan = fftw_mpi_plan_dft_2d(N0,N1,data_out,data_out2,PETSC_COMM_WORLD,FFTW_BACKWARD,FFTW_ESTIMATE); 67c4762a1bSJed Brown 689566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 699566063dSJacob Faibussowitsch if (view) PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown fftw_execute(fplan); 729566063dSJacob Faibussowitsch if (view) PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); 73c4762a1bSJed Brown 74c4762a1bSJed Brown fftw_execute(bplan); 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */ 77c4762a1bSJed Brown a = 1.0/(PetscReal)N; 789566063dSJacob Faibussowitsch PetscCall(VecScale(z,a)); 799566063dSJacob Faibussowitsch if (view) PetscCall(VecView(z, PETSC_VIEWER_STDOUT_WORLD)); 809566063dSJacob Faibussowitsch PetscCall(VecAXPY(z,-1.0,x)); 819566063dSJacob Faibussowitsch PetscCall(VecNorm(z,NORM_1,&enorm)); 82dd400576SPatrick Sanan if (enorm > 1.e-11 && rank == 0) { 839566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," Error norm of |x - z| %g\n",(double)enorm)); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* Free spaces */ 87c4762a1bSJed Brown fftw_destroy_plan(fplan); 88c4762a1bSJed Brown fftw_destroy_plan(bplan); 899566063dSJacob Faibussowitsch fftw_free(data_in); PetscCall(VecDestroy(&x)); 909566063dSJacob Faibussowitsch fftw_free(data_out); PetscCall(VecDestroy(&y)); 919566063dSJacob Faibussowitsch fftw_free(data_out2);PetscCall(VecDestroy(&z)); 92c4762a1bSJed Brown 93c4762a1bSJed Brown } else { 94c4762a1bSJed Brown /* Use PETSc-FFTW interface */ 95c4762a1bSJed Brown /*-------------------------------------------*/ 96c4762a1bSJed Brown PetscInt i,*dim,k; 97c4762a1bSJed Brown Mat A; 98c4762a1bSJed Brown 99c4762a1bSJed Brown N=1; 100c4762a1bSJed Brown for (i=1; i<5; i++) { 101c4762a1bSJed Brown DIM = i; 1029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(i,&dim)); 103c4762a1bSJed Brown for (k=0; k<i; k++) { 104c4762a1bSJed Brown dim[k]=30; 105c4762a1bSJed Brown } 106c4762a1bSJed Brown N *= dim[i-1]; 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* Create FFTW object */ 109dd400576SPatrick Sanan if (rank == 0) printf("Use PETSc-FFTW interface...%d-DIM: %d\n",(int)DIM,(int)N); 110c4762a1bSJed Brown 1119566063dSJacob Faibussowitsch PetscCall(MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A)); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* Create vectors that are compatible with parallel layout of A - must call MatCreateVecs()! */ 114c4762a1bSJed Brown 1159566063dSJacob Faibussowitsch PetscCall(MatCreateVecsFFTW(A,&x,&y,&z)); 1169566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) x, "Real space vector")); 1179566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) y, "Frequency space vector")); 1189566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) z, "Reconstructed vector")); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* Set values of space vector x */ 1219566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,rdm)); 122c4762a1bSJed Brown 1239566063dSJacob Faibussowitsch if (view) PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* Apply FFTW_FORWARD and FFTW_BACKWARD */ 1269566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,y)); 1279566063dSJacob Faibussowitsch if (view) PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); 128c4762a1bSJed Brown 1299566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A,y,z)); 130c4762a1bSJed Brown 131c4762a1bSJed Brown /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */ 132c4762a1bSJed Brown a = 1.0/(PetscReal)N; 1339566063dSJacob Faibussowitsch PetscCall(VecScale(z,a)); 1349566063dSJacob Faibussowitsch if (view) PetscCall(VecView(z,PETSC_VIEWER_STDOUT_WORLD)); 1359566063dSJacob Faibussowitsch PetscCall(VecAXPY(z,-1.0,x)); 1369566063dSJacob Faibussowitsch PetscCall(VecNorm(z,NORM_1,&enorm)); 137dd400576SPatrick Sanan if (enorm > 1.e-9 && rank == 0) { 1389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," Error norm of |x - z| %e\n",enorm)); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 1419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&z)); 1449566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 145c4762a1bSJed Brown 1469566063dSJacob Faibussowitsch PetscCall(PetscFree(dim)); 147c4762a1bSJed Brown } 148c4762a1bSJed Brown } 149c4762a1bSJed Brown 1509566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 1519566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 152b122ec5aSJacob Faibussowitsch return 0; 153c4762a1bSJed Brown } 154c4762a1bSJed Brown 155c4762a1bSJed Brown /*TEST 156c4762a1bSJed Brown 157c4762a1bSJed Brown build: 1580cf2e031SBarry Smith requires: !mpiuni fftw complex 159c4762a1bSJed Brown 160c4762a1bSJed Brown test: 161c4762a1bSJed Brown output_file: output/ex143.out 162c4762a1bSJed Brown 163c4762a1bSJed Brown test: 164c4762a1bSJed Brown suffix: 2 165c4762a1bSJed Brown nsize: 3 166c4762a1bSJed Brown output_file: output/ex143.out 167c4762a1bSJed Brown 168c4762a1bSJed Brown TEST*/ 169