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 8d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 9d71ae5a4SJacob Faibussowitsch { 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 PetscInt *indx3, tempindx, low, *indx4, tempindx1; 27c4762a1bSJed Brown 28327415f7SBarry Smith PetscFunctionBeginUser; 29c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 32c4762a1bSJed Brown 339566063dSJacob Faibussowitsch PetscCall(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 */ 589566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, n, N, (PetscScalar *)in1, &fin)); 599566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, n, N, (PetscScalar *)out, &fout)); 609566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, n, N, (PetscScalar *)in2, &fout1)); 61c4762a1bSJed Brown 62c4762a1bSJed Brown /* Set the vector with random data */ 639566063dSJacob Faibussowitsch PetscCall(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; 719566063dSJacob Faibussowitsch PetscCall(VecSetValues(fin, 1, &i, &one, INSERT_VALUES)); 72c4762a1bSJed Brown i = 1; 739566063dSJacob Faibussowitsch PetscCall(VecSetValues(fin, 1, &i, &two, INSERT_VALUES)); 74c4762a1bSJed Brown i = 4; 759566063dSJacob Faibussowitsch PetscCall(VecSetValues(fin, 1, &i, &three, INSERT_VALUES)); 76c4762a1bSJed Brown i = 5; 779566063dSJacob Faibussowitsch PetscCall(VecSetValues(fin, 1, &i, &four, INSERT_VALUES)); 789566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(fin)); 799566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(fin)); 80c4762a1bSJed Brown 819566063dSJacob Faibussowitsch PetscCall(VecSet(fout, zero)); 829566063dSJacob Faibussowitsch PetscCall(VecSet(fout1, zero)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* Get the meaningful portion of array */ 859566063dSJacob Faibussowitsch PetscCall(VecGetArray(fin, &x_arr)); 869566063dSJacob Faibussowitsch PetscCall(VecGetArray(fout1, &z_arr)); 879566063dSJacob Faibussowitsch PetscCall(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 959566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(fin, &x_arr)); 969566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(fout1, &z_arr)); 979566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(fout, &y_arr)); 98c4762a1bSJed Brown 99c4762a1bSJed Brown /* VecView(fin,PETSC_VIEWER_STDOUT_WORLD); */ 1009566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &ini)); 1019566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &final)); 1029566063dSJacob Faibussowitsch PetscCall(VecSetSizes(ini, local_n0 * N1, N0 * N1)); 1039566063dSJacob Faibussowitsch PetscCall(VecSetSizes(final, local_n0 * N1, N0 * N1)); 1049566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(ini)); 1059566063dSJacob Faibussowitsch PetscCall(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); */ 1139566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(fin, &low, NULL)); 114c4762a1bSJed Brown /*printf("The local index is %d from %d\n",low,rank); */ 1159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(local_n0 * N1, &indx3)); 1169566063dSJacob Faibussowitsch PetscCall(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 1299566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(local_n0 * N1, &x_arr, local_n0 * N1, &y_arr)); /* arr must be allocated for VecGetValues() */ 1309566063dSJacob Faibussowitsch PetscCall(VecGetValues(fin, local_n0 * N1, indx4, (PetscScalar *)x_arr)); 1319566063dSJacob Faibussowitsch PetscCall(VecSetValues(ini, local_n0 * N1, indx3, x_arr, INSERT_VALUES)); 132c4762a1bSJed Brown 1339566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(ini)); 1349566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(ini)); 135c4762a1bSJed Brown 1369566063dSJacob Faibussowitsch PetscCall(VecGetValues(fout1, local_n0 * N1, indx4, y_arr)); 1379566063dSJacob Faibussowitsch PetscCall(VecSetValues(final, local_n0 * N1, indx3, y_arr, INSERT_VALUES)); 1389566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(final)); 1399566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(final)); 1409566063dSJacob Faibussowitsch PetscCall(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; 1599566063dSJacob Faibussowitsch PetscCall(VecScale(fout1, a)); 1609566063dSJacob Faibussowitsch PetscCall(VecScale(final, a)); 161c4762a1bSJed Brown 162c4762a1bSJed Brown /* VecView(ini,PETSC_VIEWER_STDOUT_WORLD); */ 163c4762a1bSJed Brown /* VecView(final,PETSC_VIEWER_STDOUT_WORLD); */ 1649566063dSJacob Faibussowitsch PetscCall(VecAXPY(final, -1.0, ini)); 165c4762a1bSJed Brown 1669566063dSJacob Faibussowitsch PetscCall(VecNorm(final, NORM_1, &enorm)); 16748a46eb9SPierre Jolivet if (enorm > 1.e-10) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Error norm of |x - z| = %e\n", enorm)); 168c4762a1bSJed Brown 169c4762a1bSJed Brown /* Execute fftw with function fftw_execute and destroy it after execution */ 170c4762a1bSJed Brown fftw_destroy_plan(fplan); 171c4762a1bSJed Brown fftw_destroy_plan(bplan); 1729371c9d4SSatish Balay fftw_free(in1); 1739371c9d4SSatish Balay PetscCall(VecDestroy(&fin)); 1749371c9d4SSatish Balay fftw_free(out); 1759371c9d4SSatish Balay PetscCall(VecDestroy(&fout)); 1769371c9d4SSatish Balay fftw_free(in2); 1779371c9d4SSatish Balay PetscCall(VecDestroy(&fout1)); 178c4762a1bSJed Brown 1799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ini)); 1809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&final)); 181c4762a1bSJed Brown 1829566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd)); 1839566063dSJacob Faibussowitsch PetscCall(PetscFree(indx3)); 1849566063dSJacob Faibussowitsch PetscCall(PetscFree(indx4)); 1859566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 186b122ec5aSJacob Faibussowitsch return 0; 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 189c4762a1bSJed Brown /*TEST 190c4762a1bSJed Brown 191c4762a1bSJed Brown build: 1920cf2e031SBarry Smith requires: !mpiuni fftw !complex 193c4762a1bSJed Brown 194c4762a1bSJed Brown test: 195*3886731fSPierre Jolivet output_file: output/empty.out 196c4762a1bSJed Brown 197c4762a1bSJed Brown test: 198c4762a1bSJed Brown suffix: 2 199c4762a1bSJed Brown nsize: 3 200*3886731fSPierre Jolivet output_file: output/empty.out 201c4762a1bSJed Brown 202c4762a1bSJed Brown TEST*/ 203