Actual source code: ex143.c

petsc-3.3-p6 2013-02-11
  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: PetscInt main(PetscInt argc,char **args)
 19: {
 20:   PetscErrorCode  ierr;
 21:   PetscMPIInt     rank,size;
 22:   PetscInt        N0=50,N1=20,N=N0*N1;
 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, PETSC_NULL, "FFTW Options", "ex143");
 35:     PetscOptionsBool("-vec_view_draw", "View the vectors", "ex143", view, &view, PETSC_NULL);
 36:     PetscOptionsBool("-use_FFTW_interface", "Use PETSc-FFTW interface", "ex143",use_interface, &use_interface, PETSC_NULL);
 37:   PetscOptionsEnd();

 39:   PetscOptionsGetBool(PETSC_NULL,"-use_FFTW_interface",&use_interface,PETSC_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;

 53:     if (!rank) printf("Use FFTW without PETSc-FFTW interface\n");
 54:     fftw_mpi_init();
 55:     N = N0*N1;
 56:     alloc_local = fftw_mpi_local_size_2d(N0,N1,PETSC_COMM_WORLD,&local_n0,&local_0_start);

 58:     data_in   = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
 59:     data_out  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
 60:     data_out2 = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
 61:     VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_in,&x);
 62:     PetscObjectSetName((PetscObject) x, "Real Space vector");
 63:     VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out,&y);
 64:     PetscObjectSetName((PetscObject) y, "Frequency space vector");
 65:     VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out2,&z);
 66:     PetscObjectSetName((PetscObject) z, "Reconstructed vector");

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

 71:     VecSetRandom(x, rdm);
 72:     if (view){VecView(x,PETSC_VIEWER_STDOUT_WORLD);}

 74:     fftw_execute(fplan);
 75:     if (view){VecView(y,PETSC_VIEWER_STDOUT_WORLD);}

 77:     fftw_execute(bplan);

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

 89:     /* Free spaces */
 90:     fftw_destroy_plan(fplan);
 91:     fftw_destroy_plan(bplan);
 92:     fftw_free(data_in);  VecDestroy(&x);
 93:     fftw_free(data_out); VecDestroy(&y);
 94:     fftw_free(data_out2);VecDestroy(&z);

 96:   } else {
 97:     /* Use PETSc-FFTW interface                  */
 98:     /*-------------------------------------------*/
 99:     PetscInt i,*dim,k,DIM;
100:     Mat      A;

102:     N=1;
103:     for (i=1; i<5; i++){
104:       DIM = i;
105:       PetscMalloc(i*sizeof(PetscInt),&dim);
106:       for(k=0;k<i;k++){
107:         dim[k]=30;
108:       }
109:       N *= dim[i-1];

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

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

117:       /* Create vectors that are compatible with parallel layout of A - must call MatGetVecs()! */
118: 
119:       MatGetVecsFFTW(A,&x,&y,&z);
120:       PetscObjectSetName((PetscObject) x, "Real space vector");
121:       PetscObjectSetName((PetscObject) y, "Frequency space vector");
122:       PetscObjectSetName((PetscObject) z, "Reconstructed vector");

124:       /* Set values of space vector x */
125:       VecSetRandom(x,rdm);

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

129:       // Apply FFTW_FORWARD and FFTW_BACKWARD
130:       MatMult(A,x,y);
131:       if (view){VecView(y,PETSC_VIEWER_STDOUT_WORLD);}

133:       MatMultTranspose(A,y,z);

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

150:       PetscFree(dim);
151:     }
152:   }
153: 
154:   PetscRandomDestroy(&rdm);
155:   PetscFinalize();
156:   return 0;
157: }