Actual source code: ex143.c

petsc-3.10.5 2019-03-28
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>

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

 27:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 28:   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex143");
 29:   PetscOptionsBool("-vec_view draw", "View the vectors", "ex143", view, &view, NULL);
 30:   PetscOptionsBool("-use_FFTW_interface", "Use PETSc-FFTW interface", "ex143",use_interface, &use_interface, NULL);
 31:   PetscOptionsEnd();

 33:   PetscOptionsGetBool(NULL,NULL,"-use_FFTW_interface",&use_interface,NULL);
 34:   MPI_Comm_size(PETSC_COMM_WORLD, &size);
 35:   MPI_Comm_rank(PETSC_COMM_WORLD, &rank);

 37:   PetscRandomCreate(PETSC_COMM_WORLD, &rdm);
 38:   PetscRandomSetFromOptions(rdm);

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

 55:     data_in   = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
 56:     data_out  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
 57:     data_out2 = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);

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

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

 69:     VecSetRandom(x, rdm);
 70:     if (view) {VecView(x,PETSC_VIEWER_STDOUT_WORLD);}

 72:     fftw_execute(fplan);
 73:     if (view) {VecView(y,PETSC_VIEWER_STDOUT_WORLD);}

 75:     fftw_execute(bplan);

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

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

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

100:     N=1;
101:     for (i=1; i<5; i++) {
102:       DIM  = i;
103:       PetscMalloc1(i,&dim);
104:       for (k=0; k<i; k++) {
105:         dim[k]=30;
106:       }
107:       N *= dim[i-1];


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

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

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

117:       MatCreateVecsFFTW(A,&x,&y,&z);
118:       PetscObjectSetName((PetscObject) x, "Real space vector");
119:       PetscObjectSetName((PetscObject) y, "Frequency space vector");
120:       PetscObjectSetName((PetscObject) z, "Reconstructed vector");

122:       /* Set values of space vector x */
123:       VecSetRandom(x,rdm);

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

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

131:       MatMultTranspose(A,y,z);

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

143:       VecDestroy(&x);
144:       VecDestroy(&y);
145:       VecDestroy(&z);
146:       MatDestroy(&A);

148:       PetscFree(dim);
149:     }
150:   }

152:   PetscRandomDestroy(&rdm);
153:   PetscFinalize();
154:   return ierr;
155: }


158: /*TEST

160:    build:
161:       requires: fftw complex

163:    test:
164:       output_file: output/ex143.out

166:    test:
167:       suffix: 2
168:       nsize: 3
169:       output_file: output/ex143.out

171: TEST*/