xref: /petsc/src/mat/tests/ex143.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown static char help[] = "Illustrate how to use mpi FFTW and PETSc-FFTW interface \n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown   Compiling the code:
5c4762a1bSJed Brown       This code uses the complex numbers version of PETSc, so configure
6c4762a1bSJed Brown       must be run to enable this
7c4762a1bSJed Brown 
8c4762a1bSJed Brown  Usage:
9c4762a1bSJed Brown    mpiexec -n <np> ./ex143 -use_FFTW_interface NO
10c4762a1bSJed Brown    mpiexec -n <np> ./ex143 -use_FFTW_interface YES
11c4762a1bSJed Brown */
12c4762a1bSJed Brown 
13c4762a1bSJed Brown #include <petscmat.h>
14c4762a1bSJed Brown #include <fftw3-mpi.h>
15c4762a1bSJed Brown 
16c4762a1bSJed Brown int main(int argc,char **args)
17c4762a1bSJed Brown {
18c4762a1bSJed Brown   PetscMPIInt    rank,size;
19c4762a1bSJed Brown   PetscInt       N0=50,N1=20,N=N0*N1,DIM;
20c4762a1bSJed Brown   PetscRandom    rdm;
21c4762a1bSJed Brown   PetscScalar    a;
22c4762a1bSJed Brown   PetscReal      enorm;
23c4762a1bSJed Brown   Vec            x,y,z;
24c4762a1bSJed Brown   PetscBool      view=PETSC_FALSE,use_interface=PETSC_TRUE;
25c4762a1bSJed Brown 
26*327415f7SBarry Smith   PetscFunctionBeginUser;
279566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
28d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex143");
299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-vec_view draw", "View the vectors", "ex143", view, &view, NULL));
309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-use_FFTW_interface", "Use PETSc-FFTW interface", "ex143",use_interface, &use_interface, NULL));
31d0609cedSBarry Smith   PetscOptionsEnd();
32c4762a1bSJed Brown 
339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-use_FFTW_interface",&use_interface,NULL));
349566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
36c4762a1bSJed Brown 
379566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
389566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   if (!use_interface) {
41c4762a1bSJed Brown     /* Use mpi FFTW without PETSc-FFTW interface, 2D case only */
42c4762a1bSJed Brown     /*---------------------------------------------------------*/
43c4762a1bSJed Brown     fftw_plan    fplan,bplan;
44c4762a1bSJed Brown     fftw_complex *data_in,*data_out,*data_out2;
45c4762a1bSJed Brown     ptrdiff_t    alloc_local,local_n0,local_0_start;
46c4762a1bSJed Brown 
47c4762a1bSJed Brown     DIM = 2;
48dd400576SPatrick Sanan     if (rank == 0) {
499566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF,"Use FFTW without PETSc-FFTW interface, DIM %" PetscInt_FMT "\n",DIM));
50c4762a1bSJed Brown     }
51c4762a1bSJed Brown     fftw_mpi_init();
52c4762a1bSJed Brown     N           = N0*N1;
53c4762a1bSJed Brown     alloc_local = fftw_mpi_local_size_2d(N0,N1,PETSC_COMM_WORLD,&local_n0,&local_0_start);
54c4762a1bSJed Brown 
55c4762a1bSJed Brown     data_in   = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
56c4762a1bSJed Brown     data_out  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
57c4762a1bSJed Brown     data_out2 = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
58c4762a1bSJed Brown 
599566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_in,&x));
609566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) x, "Real Space vector"));
619566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out,&y));
629566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) y, "Frequency space vector"));
639566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out2,&z));
649566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) z, "Reconstructed vector"));
65c4762a1bSJed Brown 
66c4762a1bSJed Brown     fplan = fftw_mpi_plan_dft_2d(N0,N1,data_in,data_out,PETSC_COMM_WORLD,FFTW_FORWARD,FFTW_ESTIMATE);
67c4762a1bSJed Brown     bplan = fftw_mpi_plan_dft_2d(N0,N1,data_out,data_out2,PETSC_COMM_WORLD,FFTW_BACKWARD,FFTW_ESTIMATE);
68c4762a1bSJed Brown 
699566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x, rdm));
709566063dSJacob Faibussowitsch     if (view) PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
71c4762a1bSJed Brown 
72c4762a1bSJed Brown     fftw_execute(fplan);
739566063dSJacob Faibussowitsch     if (view) PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
74c4762a1bSJed Brown 
75c4762a1bSJed Brown     fftw_execute(bplan);
76c4762a1bSJed Brown 
77c4762a1bSJed Brown     /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */
78c4762a1bSJed Brown     a    = 1.0/(PetscReal)N;
799566063dSJacob Faibussowitsch     PetscCall(VecScale(z,a));
809566063dSJacob Faibussowitsch     if (view) PetscCall(VecView(z, PETSC_VIEWER_STDOUT_WORLD));
819566063dSJacob Faibussowitsch     PetscCall(VecAXPY(z,-1.0,x));
829566063dSJacob Faibussowitsch     PetscCall(VecNorm(z,NORM_1,&enorm));
83dd400576SPatrick Sanan     if (enorm > 1.e-11 && rank == 0) {
849566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %g\n",(double)enorm));
85c4762a1bSJed Brown     }
86c4762a1bSJed Brown 
87c4762a1bSJed Brown     /* Free spaces */
88c4762a1bSJed Brown     fftw_destroy_plan(fplan);
89c4762a1bSJed Brown     fftw_destroy_plan(bplan);
909566063dSJacob Faibussowitsch     fftw_free(data_in);  PetscCall(VecDestroy(&x));
919566063dSJacob Faibussowitsch     fftw_free(data_out); PetscCall(VecDestroy(&y));
929566063dSJacob Faibussowitsch     fftw_free(data_out2);PetscCall(VecDestroy(&z));
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   } else {
95c4762a1bSJed Brown     /* Use PETSc-FFTW interface                  */
96c4762a1bSJed Brown     /*-------------------------------------------*/
97c4762a1bSJed Brown     PetscInt i,*dim,k;
98c4762a1bSJed Brown     Mat      A;
99c4762a1bSJed Brown 
100c4762a1bSJed Brown     N=1;
101c4762a1bSJed Brown     for (i=1; i<5; i++) {
102c4762a1bSJed Brown       DIM  = i;
1039566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(i,&dim));
104c4762a1bSJed Brown       for (k=0; k<i; k++) {
105c4762a1bSJed Brown         dim[k]=30;
106c4762a1bSJed Brown       }
107c4762a1bSJed Brown       N *= dim[i-1];
108c4762a1bSJed Brown 
109c4762a1bSJed Brown       /* Create FFTW object */
110dd400576SPatrick Sanan       if (rank == 0) printf("Use PETSc-FFTW interface...%d-DIM: %d\n",(int)DIM,(int)N);
111c4762a1bSJed Brown 
1129566063dSJacob Faibussowitsch       PetscCall(MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A));
113c4762a1bSJed Brown 
114c4762a1bSJed Brown       /* Create vectors that are compatible with parallel layout of A - must call MatCreateVecs()! */
115c4762a1bSJed Brown 
1169566063dSJacob Faibussowitsch       PetscCall(MatCreateVecsFFTW(A,&x,&y,&z));
1179566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject) x, "Real space vector"));
1189566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject) y, "Frequency space vector"));
1199566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject) z, "Reconstructed vector"));
120c4762a1bSJed Brown 
121c4762a1bSJed Brown       /* Set values of space vector x */
1229566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(x,rdm));
123c4762a1bSJed Brown 
1249566063dSJacob Faibussowitsch       if (view) PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown       /* Apply FFTW_FORWARD and FFTW_BACKWARD */
1279566063dSJacob Faibussowitsch       PetscCall(MatMult(A,x,y));
1289566063dSJacob Faibussowitsch       if (view) PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
129c4762a1bSJed Brown 
1309566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(A,y,z));
131c4762a1bSJed Brown 
132c4762a1bSJed Brown       /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */
133c4762a1bSJed Brown       a    = 1.0/(PetscReal)N;
1349566063dSJacob Faibussowitsch       PetscCall(VecScale(z,a));
1359566063dSJacob Faibussowitsch       if (view) PetscCall(VecView(z,PETSC_VIEWER_STDOUT_WORLD));
1369566063dSJacob Faibussowitsch       PetscCall(VecAXPY(z,-1.0,x));
1379566063dSJacob Faibussowitsch       PetscCall(VecNorm(z,NORM_1,&enorm));
138dd400576SPatrick Sanan       if (enorm > 1.e-9 && rank == 0) {
1399566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %e\n",enorm));
140c4762a1bSJed Brown       }
141c4762a1bSJed Brown 
1429566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&x));
1439566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&y));
1449566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&z));
1459566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&A));
146c4762a1bSJed Brown 
1479566063dSJacob Faibussowitsch       PetscCall(PetscFree(dim));
148c4762a1bSJed Brown     }
149c4762a1bSJed Brown   }
150c4762a1bSJed Brown 
1519566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
1529566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
153b122ec5aSJacob Faibussowitsch   return 0;
154c4762a1bSJed Brown }
155c4762a1bSJed Brown 
156c4762a1bSJed Brown /*TEST
157c4762a1bSJed Brown 
158c4762a1bSJed Brown    build:
1590cf2e031SBarry Smith       requires: !mpiuni fftw complex
160c4762a1bSJed Brown 
161c4762a1bSJed Brown    test:
162c4762a1bSJed Brown       output_file: output/ex143.out
163c4762a1bSJed Brown 
164c4762a1bSJed Brown    test:
165c4762a1bSJed Brown       suffix: 2
166c4762a1bSJed Brown       nsize: 3
167c4762a1bSJed Brown       output_file: output/ex143.out
168c4762a1bSJed Brown 
169c4762a1bSJed Brown TEST*/
170