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 169371c9d4SSatish Balay int main(int argc, char **args) { 17c4762a1bSJed Brown PetscMPIInt rank, size; 18c4762a1bSJed Brown PetscInt N0 = 50, N1 = 20, N = N0 * N1, DIM; 19c4762a1bSJed Brown PetscRandom rdm; 20c4762a1bSJed Brown PetscScalar a; 21c4762a1bSJed Brown PetscReal enorm; 22c4762a1bSJed Brown Vec x, y, z; 23c4762a1bSJed Brown PetscBool view = PETSC_FALSE, use_interface = PETSC_TRUE; 24c4762a1bSJed Brown 25327415f7SBarry Smith PetscFunctionBeginUser; 269566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 27d0609cedSBarry 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)); 30d0609cedSBarry 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; 47*48a46eb9SPierre Jolivet if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Use FFTW without PETSc-FFTW interface, DIM %" PetscInt_FMT "\n", DIM)); 48c4762a1bSJed Brown fftw_mpi_init(); 49c4762a1bSJed Brown N = N0 * N1; 50c4762a1bSJed Brown alloc_local = fftw_mpi_local_size_2d(N0, N1, PETSC_COMM_WORLD, &local_n0, &local_0_start); 51c4762a1bSJed Brown 52c4762a1bSJed Brown data_in = (fftw_complex *)fftw_malloc(sizeof(fftw_complex) * alloc_local); 53c4762a1bSJed Brown data_out = (fftw_complex *)fftw_malloc(sizeof(fftw_complex) * alloc_local); 54c4762a1bSJed Brown data_out2 = (fftw_complex *)fftw_malloc(sizeof(fftw_complex) * alloc_local); 55c4762a1bSJed Brown 569566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, (PetscInt)local_n0 * N1, (PetscInt)N, (const PetscScalar *)data_in, &x)); 579566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x, "Real Space vector")); 589566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, (PetscInt)local_n0 * N1, (PetscInt)N, (const PetscScalar *)data_out, &y)); 599566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)y, "Frequency space vector")); 609566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, (PetscInt)local_n0 * N1, (PetscInt)N, (const PetscScalar *)data_out2, &z)); 619566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)z, "Reconstructed vector")); 62c4762a1bSJed Brown 63c4762a1bSJed Brown fplan = fftw_mpi_plan_dft_2d(N0, N1, data_in, data_out, PETSC_COMM_WORLD, FFTW_FORWARD, FFTW_ESTIMATE); 64c4762a1bSJed Brown bplan = fftw_mpi_plan_dft_2d(N0, N1, data_out, data_out2, PETSC_COMM_WORLD, FFTW_BACKWARD, FFTW_ESTIMATE); 65c4762a1bSJed Brown 669566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 679566063dSJacob Faibussowitsch if (view) PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown fftw_execute(fplan); 709566063dSJacob Faibussowitsch if (view) PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD)); 71c4762a1bSJed Brown 72c4762a1bSJed Brown fftw_execute(bplan); 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */ 75c4762a1bSJed Brown a = 1.0 / (PetscReal)N; 769566063dSJacob Faibussowitsch PetscCall(VecScale(z, a)); 779566063dSJacob Faibussowitsch if (view) PetscCall(VecView(z, PETSC_VIEWER_STDOUT_WORLD)); 789566063dSJacob Faibussowitsch PetscCall(VecAXPY(z, -1.0, x)); 799566063dSJacob Faibussowitsch PetscCall(VecNorm(z, NORM_1, &enorm)); 80*48a46eb9SPierre Jolivet if (enorm > 1.e-11 && rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Error norm of |x - z| %g\n", (double)enorm)); 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* Free spaces */ 83c4762a1bSJed Brown fftw_destroy_plan(fplan); 84c4762a1bSJed Brown fftw_destroy_plan(bplan); 859371c9d4SSatish Balay fftw_free(data_in); 869371c9d4SSatish Balay PetscCall(VecDestroy(&x)); 879371c9d4SSatish Balay fftw_free(data_out); 889371c9d4SSatish Balay PetscCall(VecDestroy(&y)); 899371c9d4SSatish Balay fftw_free(data_out2); 909371c9d4SSatish Balay PetscCall(VecDestroy(&z)); 91c4762a1bSJed Brown 92c4762a1bSJed Brown } else { 93c4762a1bSJed Brown /* Use PETSc-FFTW interface */ 94c4762a1bSJed Brown /*-------------------------------------------*/ 95c4762a1bSJed Brown PetscInt i, *dim, k; 96c4762a1bSJed Brown Mat A; 97c4762a1bSJed Brown 98c4762a1bSJed Brown N = 1; 99c4762a1bSJed Brown for (i = 1; i < 5; i++) { 100c4762a1bSJed Brown DIM = i; 1019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(i, &dim)); 1029371c9d4SSatish Balay for (k = 0; k < i; k++) { dim[k] = 30; } 103c4762a1bSJed Brown N *= dim[i - 1]; 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* Create FFTW object */ 106dd400576SPatrick Sanan if (rank == 0) printf("Use PETSc-FFTW interface...%d-DIM: %d\n", (int)DIM, (int)N); 107c4762a1bSJed Brown 1089566063dSJacob Faibussowitsch PetscCall(MatCreateFFT(PETSC_COMM_WORLD, DIM, dim, MATFFTW, &A)); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* Create vectors that are compatible with parallel layout of A - must call MatCreateVecs()! */ 111c4762a1bSJed Brown 1129566063dSJacob Faibussowitsch PetscCall(MatCreateVecsFFTW(A, &x, &y, &z)); 1139566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x, "Real space vector")); 1149566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)y, "Frequency space vector")); 1159566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)z, "Reconstructed vector")); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* Set values of space vector x */ 1189566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 119c4762a1bSJed Brown 1209566063dSJacob Faibussowitsch if (view) PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* Apply FFTW_FORWARD and FFTW_BACKWARD */ 1239566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, y)); 1249566063dSJacob Faibussowitsch if (view) PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD)); 125c4762a1bSJed Brown 1269566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, y, z)); 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */ 129c4762a1bSJed Brown a = 1.0 / (PetscReal)N; 1309566063dSJacob Faibussowitsch PetscCall(VecScale(z, a)); 1319566063dSJacob Faibussowitsch if (view) PetscCall(VecView(z, PETSC_VIEWER_STDOUT_WORLD)); 1329566063dSJacob Faibussowitsch PetscCall(VecAXPY(z, -1.0, x)); 1339566063dSJacob Faibussowitsch PetscCall(VecNorm(z, NORM_1, &enorm)); 134*48a46eb9SPierre Jolivet if (enorm > 1.e-9 && rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Error norm of |x - z| %e\n", enorm)); 135c4762a1bSJed Brown 1369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&z)); 1399566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 140c4762a1bSJed Brown 1419566063dSJacob Faibussowitsch PetscCall(PetscFree(dim)); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 1459566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 1469566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 147b122ec5aSJacob Faibussowitsch return 0; 148c4762a1bSJed Brown } 149c4762a1bSJed Brown 150c4762a1bSJed Brown /*TEST 151c4762a1bSJed Brown 152c4762a1bSJed Brown build: 1530cf2e031SBarry Smith requires: !mpiuni fftw complex 154c4762a1bSJed Brown 155c4762a1bSJed Brown test: 156c4762a1bSJed Brown output_file: output/ex143.out 157c4762a1bSJed Brown 158c4762a1bSJed Brown test: 159c4762a1bSJed Brown suffix: 2 160c4762a1bSJed Brown nsize: 3 161c4762a1bSJed Brown output_file: output/ex143.out 162c4762a1bSJed Brown 163c4762a1bSJed Brown TEST*/ 164