Actual source code: ex158.c
petsc-3.7.3 2016-08-01
1: static char help[] = "Illustrate how to use mpi FFTW and PETSc-FFTW interface \n\n";
3: /*
4: Usage:
5: mpiexec -n <np> ./ex158 -use_FFTW_interface NO
6: mpiexec -n <np> ./ex158 -use_FFTW_interface YES
7: */
9: #include <petscmat.h>
10: #include <fftw3-mpi.h>
14: int main(int argc,char **args)
15: {
17: PetscMPIInt rank,size;
18: PetscInt N0=50,N1=20,N=N0*N1;
19: PetscRandom rdm;
20: PetscScalar a;
21: PetscReal enorm;
22: Vec x,y,z;
23: PetscBool view=PETSC_FALSE,use_interface=PETSC_TRUE;
25: PetscInitialize(&argc,&args,(char*)0,help);
26: #if defined(PETSC_USE_COMPLEX)
27: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers. Your current scalar type is complex");
28: #endif
30: PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex158");
31: PetscOptionsBool("-use_FFTW_interface", "Use PETSc-FFTW interface", "ex158",use_interface, &use_interface, NULL);
32: PetscOptionsEnd();
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;
47: if (!rank) printf("Use FFTW without PETSc-FFTW interface\n");
48: fftw_mpi_init();
49: N = N0*N1;
50: alloc_local = fftw_mpi_local_size_2d(N0,N1,PETSC_COMM_WORLD,&local_n0,&local_0_start);
52: data_in = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
53: data_out = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
54: data_out2 = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
56: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_in,&x);
57: PetscObjectSetName((PetscObject) x, "Real Space vector");
58: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out,&y);
59: PetscObjectSetName((PetscObject) y, "Frequency space vector");
60: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out2,&z);
61: PetscObjectSetName((PetscObject) z, "Reconstructed vector");
63: fplan = fftw_mpi_plan_dft_2d(N0,N1,data_in,data_out,PETSC_COMM_WORLD,FFTW_FORWARD,FFTW_ESTIMATE);
64: bplan = fftw_mpi_plan_dft_2d(N0,N1,data_out,data_out2,PETSC_COMM_WORLD,FFTW_BACKWARD,FFTW_ESTIMATE);
66: VecSetRandom(x, rdm);
67: if (view) {VecView(x,PETSC_VIEWER_STDOUT_WORLD);}
69: fftw_execute(fplan);
70: if (view) {VecView(y,PETSC_VIEWER_STDOUT_WORLD);}
72: fftw_execute(bplan);
74: /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */
75: a = 1.0/(PetscReal)N;
76: VecScale(z,a);
77: if (view) {VecView(z, PETSC_VIEWER_STDOUT_WORLD);}
78: VecAXPY(z,-1.0,x);
79: VecNorm(z,NORM_1,&enorm);
80: if (enorm > 1.e-11) {
81: PetscPrintf(PETSC_COMM_SELF," Error norm of |x - z| %g\n",(double)enorm);
82: }
84: /* Free spaces */
85: fftw_destroy_plan(fplan);
86: fftw_destroy_plan(bplan);
87: fftw_free(data_in); VecDestroy(&x);
88: fftw_free(data_out); VecDestroy(&y);
89: fftw_free(data_out2);VecDestroy(&z);
91: } else {
92: /* Use PETSc-FFTW interface */
93: /*-------------------------------------------*/
94: PetscInt i,*dim,k,DIM;
95: Mat A;
96: Vec input,output;
98: N=30;
99: for (i=2; i<3; i++) { /* (i=3,4: -- error in VecScatterPetscToFFTW(A,input,x); */
100: DIM = i;
101: PetscMalloc1(i,&dim);
102: for (k=0; k<i; k++) {
103: dim[k]=30;
104: }
105: N *= dim[i-1];
107: /* Create FFTW object */
108: if (!rank) {
109: PetscPrintf(PETSC_COMM_SELF,"Use PETSc-FFTW interface...%d-DIM:%d \n",DIM,N);
110: }
111: MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A);
113: /* Create FFTW vectors that are compatible with parallel layout of A */
114: MatCreateVecsFFTW(A,&x,&y,&z);
115: PetscObjectSetName((PetscObject) x, "Real space vector");
116: PetscObjectSetName((PetscObject) y, "Frequency space vector");
117: PetscObjectSetName((PetscObject) z, "Reconstructed vector");
119: /* Create and set PETSc vector */
120: VecCreate(PETSC_COMM_WORLD,&input);
121: VecSetSizes(input,PETSC_DECIDE,N);
122: VecSetFromOptions(input);
123: VecSetRandom(input,rdm);
124: VecDuplicate(input,&output);
125: if (view) {VecView(input,PETSC_VIEWER_STDOUT_WORLD);}
127: /* Vector input is copied to another vector x using VecScatterPetscToFFTW. This is because the user data
128: can have any parallel layout. But FFTW requires special parallel layout of the data. Hence the original
129: data which is in the vector "input" here, needs to be copied to a vector x, which has the correct parallel
130: layout for FFTW. Also, during parallel real transform, this pads extra zeros automatically
131: at the end of last dimension. This padding is required by FFTW to perform parallel real D.F.T. */
132: VecScatterPetscToFFTW(A,input,x);/* buggy for dim = 3, 4... */
134: /* Apply FFTW_FORWARD and FFTW_BACKWARD */
135: MatMult(A,x,y);
136: if (view) {VecView(y,PETSC_VIEWER_STDOUT_WORLD);}
137: MatMultTranspose(A,y,z);
139: /* Output from Backward DFT needs to be modified to obtain user readable data the routine VecScatterFFTWToPetsc
140: performs the job. In some sense this is the reverse operation of VecScatterPetscToFFTW. This routine gets rid of
141: the extra spaces that were artificially padded to perform real parallel transform. */
142: VecScatterFFTWToPetsc(A,z,output);
144: /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */
145: a = 1.0/(PetscReal)N;
146: VecScale(output,a);
147: if (view) {VecView(output,PETSC_VIEWER_STDOUT_WORLD);}
148: VecAXPY(output,-1.0,input);
149: VecNorm(output,NORM_1,&enorm);
150: if (enorm > 1.e-09 && !rank) {
151: PetscPrintf(PETSC_COMM_SELF," Error norm of |x - z| %e\n",enorm);
152: }
154: /* Free spaces */
155: PetscFree(dim);
156: VecDestroy(&input);
157: VecDestroy(&output);
158: VecDestroy(&x);
159: VecDestroy(&y);
160: VecDestroy(&z);
161: MatDestroy(&A);
162: }
163: }
164: PetscRandomDestroy(&rdm);
165: PetscFinalize();
166: return 0;
167: }