Actual source code: ex143.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  1: static char help[] = "Illustrate how to use mpi FFTW and PETSc-FFTW interface \n\n";

  3: /*
  4:   Compiling the code:
  5:       This code uses the complex numbers version of PETSc, so configure
  6:       must be run to enable this

  8:  Usage:
  9:    mpiexec -n <np> ./ex143 -use_FFTW_interface NO
 10:    mpiexec -n <np> ./ex143 -use_FFTW_interface YES
 11: */

 13: #include <petscmat.h>
 14: #include <fftw3-mpi.h>

 18: int main(int argc,char **args)
 19: {
 21:   PetscMPIInt    rank,size;
 22:   PetscInt       N0=50,N1=20,N=N0*N1,DIM;
 23:   PetscRandom    rdm;
 24:   PetscScalar    a;
 25:   PetscReal      enorm;
 26:   Vec            x,y,z;
 27:   PetscBool      view=PETSC_FALSE,use_interface=PETSC_TRUE;

 29:   PetscInitialize(&argc,&args,(char*)0,help);
 30: #if !defined(PETSC_USE_COMPLEX)
 31:   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires complex numbers");
 32: #endif

 34:   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex143");
 35:   PetscOptionsBool("-vec_view draw", "View the vectors", "ex143", view, &view, NULL);
 36:   PetscOptionsBool("-use_FFTW_interface", "Use PETSc-FFTW interface", "ex143",use_interface, &use_interface, NULL);
 37:   PetscOptionsEnd();

 39:   PetscOptionsGetBool(NULL,"-use_FFTW_interface",&use_interface,NULL);
 40:   MPI_Comm_size(PETSC_COMM_WORLD, &size);
 41:   MPI_Comm_rank(PETSC_COMM_WORLD, &rank);

 43:   PetscRandomCreate(PETSC_COMM_WORLD, &rdm);
 44:   PetscRandomSetFromOptions(rdm);

 46:   if (!use_interface) {
 47:     /* Use mpi FFTW without PETSc-FFTW interface, 2D case only */
 48:     /*---------------------------------------------------------*/
 49:     fftw_plan    fplan,bplan;
 50:     fftw_complex *data_in,*data_out,*data_out2;
 51:     ptrdiff_t    alloc_local,local_n0,local_0_start;
 52: 
 53:     DIM = 2;
 54:     if (!rank) {
 55:       PetscPrintf(PETSC_COMM_SELF,"Use FFTW without PETSc-FFTW interface, DIM %D\n",DIM);
 56:     }
 57:     fftw_mpi_init();
 58:     N           = N0*N1;
 59:     alloc_local = fftw_mpi_local_size_2d(N0,N1,PETSC_COMM_WORLD,&local_n0,&local_0_start);

 61:     data_in   = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
 62:     data_out  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
 63:     data_out2 = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);

 65:     VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_in,&x);
 66:     PetscObjectSetName((PetscObject) x, "Real Space vector");
 67:     VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out,&y);
 68:     PetscObjectSetName((PetscObject) y, "Frequency space vector");
 69:     VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out2,&z);
 70:     PetscObjectSetName((PetscObject) z, "Reconstructed vector");

 72:     fplan = fftw_mpi_plan_dft_2d(N0,N1,data_in,data_out,PETSC_COMM_WORLD,FFTW_FORWARD,FFTW_ESTIMATE);
 73:     bplan = fftw_mpi_plan_dft_2d(N0,N1,data_out,data_out2,PETSC_COMM_WORLD,FFTW_BACKWARD,FFTW_ESTIMATE);

 75:     VecSetRandom(x, rdm);
 76:     if (view) {VecView(x,PETSC_VIEWER_STDOUT_WORLD);}

 78:     fftw_execute(fplan);
 79:     if (view) {VecView(y,PETSC_VIEWER_STDOUT_WORLD);}

 81:     fftw_execute(bplan);

 83:     /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */
 84:     a    = 1.0/(PetscReal)N;
 85:     VecScale(z,a);
 86:     if (view) {VecView(z, PETSC_VIEWER_STDOUT_WORLD);}
 87:     VecAXPY(z,-1.0,x);
 88:     VecNorm(z,NORM_1,&enorm);
 89:     if (enorm > 1.e-11 && !rank) {
 90:       PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %g\n",(double)enorm);
 91:     }

 93:     /* Free spaces */
 94:     fftw_destroy_plan(fplan);
 95:     fftw_destroy_plan(bplan);
 96:     fftw_free(data_in);  VecDestroy(&x);
 97:     fftw_free(data_out); VecDestroy(&y);
 98:     fftw_free(data_out2);VecDestroy(&z);

100:   } else {
101:     /* Use PETSc-FFTW interface                  */
102:     /*-------------------------------------------*/
103:     PetscInt i,*dim,k;
104:     Mat      A;

106:     N=1;
107:     for (i=1; i<5; i++) {
108:       DIM  = i;
109:       PetscMalloc1(i,&dim);
110:       for (k=0; k<i; k++) {
111:         dim[k]=30;
112:       }
113:       N *= dim[i-1];


116:       /* Create FFTW object */
117:       if (!rank) printf("Use PETSc-FFTW interface...%d-DIM: %d\n",(int)DIM,(int)N);

119:       MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A);

121:       /* Create vectors that are compatible with parallel layout of A - must call MatCreateVecs()! */

123:       MatCreateVecsFFTW(A,&x,&y,&z);
124:       PetscObjectSetName((PetscObject) x, "Real space vector");
125:       PetscObjectSetName((PetscObject) y, "Frequency space vector");
126:       PetscObjectSetName((PetscObject) z, "Reconstructed vector");

128:       /* Set values of space vector x */
129:       VecSetRandom(x,rdm);

131:       if (view) {VecView(x,PETSC_VIEWER_STDOUT_WORLD);}

133:       /* Apply FFTW_FORWARD and FFTW_BACKWARD */
134:       MatMult(A,x,y);
135:       if (view) {VecView(y,PETSC_VIEWER_STDOUT_WORLD);}

137:       MatMultTranspose(A,y,z);

139:       /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */
140:       a    = 1.0/(PetscReal)N;
141:       VecScale(z,a);
142:       if (view) {VecView(z,PETSC_VIEWER_STDOUT_WORLD);}
143:       VecAXPY(z,-1.0,x);
144:       VecNorm(z,NORM_1,&enorm);
145:       if (enorm > 1.e-9 && !rank) {
146:         PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %e\n",enorm);
147:       }

149:       VecDestroy(&x);
150:       VecDestroy(&y);
151:       VecDestroy(&z);
152:       MatDestroy(&A);

154:       PetscFree(dim);
155:     }
156:   }

158:   PetscRandomDestroy(&rdm);
159:   PetscFinalize();
160:   return 0;
161: }