xref: /petsc/src/mat/tests/ex158.c (revision d0609ced746bc51b019815ca91d747429db24893)
1c4762a1bSJed Brown static char help[] = "Illustrate how to use mpi FFTW and PETSc-FFTW interface \n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown  Usage:
5c4762a1bSJed Brown    mpiexec -n <np> ./ex158 -use_FFTW_interface NO
6c4762a1bSJed Brown    mpiexec -n <np> ./ex158 -use_FFTW_interface YES
7c4762a1bSJed Brown */
8c4762a1bSJed Brown 
9c4762a1bSJed Brown #include <petscmat.h>
10c4762a1bSJed Brown #include <fftw3-mpi.h>
11c4762a1bSJed Brown 
12c4762a1bSJed Brown int main(int argc,char **args)
13c4762a1bSJed Brown {
14c4762a1bSJed Brown   PetscMPIInt    rank,size;
15c4762a1bSJed Brown   PetscInt       N0=50,N1=20,N=N0*N1;
16c4762a1bSJed Brown   PetscRandom    rdm;
17c4762a1bSJed Brown   PetscScalar    a;
18c4762a1bSJed Brown   PetscReal      enorm;
19c4762a1bSJed Brown   Vec            x,y,z;
20c4762a1bSJed Brown   PetscBool      view=PETSC_FALSE,use_interface=PETSC_TRUE;
21c4762a1bSJed Brown 
229566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
23c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
24c4762a1bSJed Brown   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers. Your current scalar type is complex");
25c4762a1bSJed Brown #endif
26c4762a1bSJed Brown 
27*d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex158");
289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-use_FFTW_interface", "Use PETSc-FFTW interface", "ex158",use_interface, &use_interface, NULL));
29*d0609cedSBarry Smith   PetscOptionsEnd();
30c4762a1bSJed Brown 
319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
33c4762a1bSJed Brown 
349566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
359566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   if (!use_interface) {
38c4762a1bSJed Brown     /* Use mpi FFTW without PETSc-FFTW interface, 2D case only */
39c4762a1bSJed Brown     /*---------------------------------------------------------*/
40c4762a1bSJed Brown     fftw_plan    fplan,bplan;
41c4762a1bSJed Brown     fftw_complex *data_in,*data_out,*data_out2;
42c4762a1bSJed Brown     ptrdiff_t    alloc_local,local_n0,local_0_start;
43c4762a1bSJed Brown 
44dd400576SPatrick Sanan     if (rank == 0) printf("Use FFTW without PETSc-FFTW interface\n");
45c4762a1bSJed Brown     fftw_mpi_init();
46c4762a1bSJed Brown     N           = N0*N1;
47c4762a1bSJed Brown     alloc_local = fftw_mpi_local_size_2d(N0,N1,PETSC_COMM_WORLD,&local_n0,&local_0_start);
48c4762a1bSJed Brown 
49c4762a1bSJed Brown     data_in   = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
50c4762a1bSJed Brown     data_out  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
51c4762a1bSJed Brown     data_out2 = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
52c4762a1bSJed Brown 
539566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_in,&x));
549566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) x, "Real Space vector"));
559566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out,&y));
569566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) y, "Frequency space vector"));
579566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out2,&z));
589566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) z, "Reconstructed vector"));
59c4762a1bSJed Brown 
60c4762a1bSJed Brown     fplan = fftw_mpi_plan_dft_2d(N0,N1,data_in,data_out,PETSC_COMM_WORLD,FFTW_FORWARD,FFTW_ESTIMATE);
61c4762a1bSJed Brown     bplan = fftw_mpi_plan_dft_2d(N0,N1,data_out,data_out2,PETSC_COMM_WORLD,FFTW_BACKWARD,FFTW_ESTIMATE);
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x, rdm));
649566063dSJacob Faibussowitsch     if (view) PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
65c4762a1bSJed Brown 
66c4762a1bSJed Brown     fftw_execute(fplan);
679566063dSJacob Faibussowitsch     if (view) PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
68c4762a1bSJed Brown 
69c4762a1bSJed Brown     fftw_execute(bplan);
70c4762a1bSJed Brown 
71c4762a1bSJed Brown     /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */
72c4762a1bSJed Brown     a    = 1.0/(PetscReal)N;
739566063dSJacob Faibussowitsch     PetscCall(VecScale(z,a));
749566063dSJacob Faibussowitsch     if (view) PetscCall(VecView(z, PETSC_VIEWER_STDOUT_WORLD));
759566063dSJacob Faibussowitsch     PetscCall(VecAXPY(z,-1.0,x));
769566063dSJacob Faibussowitsch     PetscCall(VecNorm(z,NORM_1,&enorm));
77c4762a1bSJed Brown     if (enorm > 1.e-11) {
789566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %g\n",(double)enorm));
79c4762a1bSJed Brown     }
80c4762a1bSJed Brown 
81c4762a1bSJed Brown     /* Free spaces */
82c4762a1bSJed Brown     fftw_destroy_plan(fplan);
83c4762a1bSJed Brown     fftw_destroy_plan(bplan);
849566063dSJacob Faibussowitsch     fftw_free(data_in);  PetscCall(VecDestroy(&x));
859566063dSJacob Faibussowitsch     fftw_free(data_out); PetscCall(VecDestroy(&y));
869566063dSJacob Faibussowitsch     fftw_free(data_out2);PetscCall(VecDestroy(&z));
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   } else {
89c4762a1bSJed Brown     /* Use PETSc-FFTW interface                  */
90c4762a1bSJed Brown     /*-------------------------------------------*/
91c4762a1bSJed Brown     PetscInt i,*dim,k,DIM;
92c4762a1bSJed Brown     Mat      A;
93c4762a1bSJed Brown     Vec      input,output;
94c4762a1bSJed Brown 
95c4762a1bSJed Brown     N=30;
96c4762a1bSJed Brown     for (i=2; i<3; i++) { /* (i=3,4: -- error in VecScatterPetscToFFTW(A,input,x); */
97c4762a1bSJed Brown       DIM  = i;
989566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(i,&dim));
99c4762a1bSJed Brown       for (k=0; k<i; k++) {
100c4762a1bSJed Brown         dim[k]=30;
101c4762a1bSJed Brown       }
102c4762a1bSJed Brown       N *= dim[i-1];
103c4762a1bSJed Brown 
104c4762a1bSJed Brown       /* Create FFTW object */
105dd400576SPatrick Sanan       if (rank == 0) {
1069566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF,"Use PETSc-FFTW interface...%d-DIM:%d \n",DIM,N));
107c4762a1bSJed Brown       }
1089566063dSJacob Faibussowitsch       PetscCall(MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A));
109c4762a1bSJed Brown 
110c4762a1bSJed Brown       /* Create FFTW vectors that are compatible with parallel layout of A */
1119566063dSJacob Faibussowitsch       PetscCall(MatCreateVecsFFTW(A,&x,&y,&z));
1129566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject) x, "Real space vector"));
1139566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject) y, "Frequency space vector"));
1149566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject) z, "Reconstructed vector"));
115c4762a1bSJed Brown 
116c4762a1bSJed Brown       /* Create and set PETSc vector */
1179566063dSJacob Faibussowitsch       PetscCall(VecCreate(PETSC_COMM_WORLD,&input));
1189566063dSJacob Faibussowitsch       PetscCall(VecSetSizes(input,PETSC_DECIDE,N));
1199566063dSJacob Faibussowitsch       PetscCall(VecSetFromOptions(input));
1209566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(input,rdm));
1219566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(input,&output));
1229566063dSJacob Faibussowitsch       if (view) PetscCall(VecView(input,PETSC_VIEWER_STDOUT_WORLD));
123c4762a1bSJed Brown 
124c4762a1bSJed Brown       /* Vector input is copied to another vector x using VecScatterPetscToFFTW. This is because the user data
125c4762a1bSJed Brown          can have any parallel layout. But FFTW requires special parallel layout of the data. Hence the original
126c4762a1bSJed Brown          data which is in the vector "input" here, needs to be copied to a vector x, which has the correct parallel
127c4762a1bSJed Brown          layout for FFTW. Also, during parallel real transform, this pads extra zeros automatically
128c4762a1bSJed Brown          at the end of last  dimension. This padding is required by FFTW to perform parallel real D.F.T.  */
1299566063dSJacob Faibussowitsch       PetscCall(VecScatterPetscToFFTW(A,input,x));/* buggy for dim = 3, 4... */
130c4762a1bSJed Brown 
131c4762a1bSJed Brown       /* Apply FFTW_FORWARD and FFTW_BACKWARD */
1329566063dSJacob Faibussowitsch       PetscCall(MatMult(A,x,y));
1339566063dSJacob Faibussowitsch       if (view) PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
1349566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(A,y,z));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown       /* Output from Backward DFT needs to be modified to obtain user readable data the routine VecScatterFFTWToPetsc
137c4762a1bSJed Brown          performs the job. In some sense this is the reverse operation of VecScatterPetscToFFTW. This routine gets rid of
138c4762a1bSJed Brown          the extra spaces that were artificially padded to perform real parallel transform.    */
1399566063dSJacob Faibussowitsch       PetscCall(VecScatterFFTWToPetsc(A,z,output));
140c4762a1bSJed Brown 
141c4762a1bSJed Brown       /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */
142c4762a1bSJed Brown       a    = 1.0/(PetscReal)N;
1439566063dSJacob Faibussowitsch       PetscCall(VecScale(output,a));
1449566063dSJacob Faibussowitsch       if (view) PetscCall(VecView(output,PETSC_VIEWER_STDOUT_WORLD));
1459566063dSJacob Faibussowitsch       PetscCall(VecAXPY(output,-1.0,input));
1469566063dSJacob Faibussowitsch       PetscCall(VecNorm(output,NORM_1,&enorm));
147dd400576SPatrick Sanan       if (enorm > 1.e-09 && rank == 0) {
1489566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %e\n",enorm));
149c4762a1bSJed Brown       }
150c4762a1bSJed Brown 
151c4762a1bSJed Brown       /* Free spaces */
1529566063dSJacob Faibussowitsch       PetscCall(PetscFree(dim));
1539566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&input));
1549566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&output));
1559566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&x));
1569566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&y));
1579566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&z));
1589566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&A));
159c4762a1bSJed Brown     }
160c4762a1bSJed Brown   }
1619566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
1629566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
163b122ec5aSJacob Faibussowitsch   return 0;
164c4762a1bSJed Brown }
165c4762a1bSJed Brown 
166c4762a1bSJed Brown /*TEST
167c4762a1bSJed Brown 
168c4762a1bSJed Brown    build:
1690cf2e031SBarry Smith       requires: !mpiuni fftw !complex
170c4762a1bSJed Brown 
171c4762a1bSJed Brown    test:
172c4762a1bSJed Brown       output_file: output/ex158.out
173c4762a1bSJed Brown 
174c4762a1bSJed Brown    test:
175c4762a1bSJed Brown       suffix: 2
176c4762a1bSJed Brown       nsize: 3
177c4762a1bSJed Brown 
178c4762a1bSJed Brown TEST*/
179