xref: /petsc/src/mat/tests/ex143.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscErrorCode ierr;
19c4762a1bSJed Brown   PetscMPIInt    rank,size;
20c4762a1bSJed Brown   PetscInt       N0=50,N1=20,N=N0*N1,DIM;
21c4762a1bSJed Brown   PetscRandom    rdm;
22c4762a1bSJed Brown   PetscScalar    a;
23c4762a1bSJed Brown   PetscReal      enorm;
24c4762a1bSJed Brown   Vec            x,y,z;
25c4762a1bSJed Brown   PetscBool      view=PETSC_FALSE,use_interface=PETSC_TRUE;
26c4762a1bSJed Brown 
27*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
28c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex143");CHKERRQ(ierr);
295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-vec_view draw", "View the vectors", "ex143", view, &view, NULL));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-use_FFTW_interface", "Use PETSc-FFTW interface", "ex143",use_interface, &use_interface, NULL));
31c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
32c4762a1bSJed Brown 
335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-use_FFTW_interface",&use_interface,NULL));
345f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
355f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
36c4762a1bSJed Brown 
375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(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) {
495f80ce2aSJacob Faibussowitsch       CHKERRQ(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 
595f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_in,&x));
605f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) x, "Real Space vector"));
615f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out,&y));
625f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) y, "Frequency space vector"));
635f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out2,&z));
645f80ce2aSJacob Faibussowitsch     CHKERRQ(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 
695f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(x, rdm));
705f80ce2aSJacob Faibussowitsch     if (view) CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
71c4762a1bSJed Brown 
72c4762a1bSJed Brown     fftw_execute(fplan);
735f80ce2aSJacob Faibussowitsch     if (view) CHKERRQ(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;
795f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScale(z,a));
805f80ce2aSJacob Faibussowitsch     if (view) CHKERRQ(VecView(z, PETSC_VIEWER_STDOUT_WORLD));
815f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(z,-1.0,x));
825f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(z,NORM_1,&enorm));
83dd400576SPatrick Sanan     if (enorm > 1.e-11 && rank == 0) {
845f80ce2aSJacob Faibussowitsch       CHKERRQ(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);
905f80ce2aSJacob Faibussowitsch     fftw_free(data_in);  CHKERRQ(VecDestroy(&x));
915f80ce2aSJacob Faibussowitsch     fftw_free(data_out); CHKERRQ(VecDestroy(&y));
925f80ce2aSJacob Faibussowitsch     fftw_free(data_out2);CHKERRQ(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;
1035f80ce2aSJacob Faibussowitsch       CHKERRQ(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 
1125f80ce2aSJacob Faibussowitsch       CHKERRQ(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 
1165f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateVecsFFTW(A,&x,&y,&z));
1175f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectSetName((PetscObject) x, "Real space vector"));
1185f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectSetName((PetscObject) y, "Frequency space vector"));
1195f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectSetName((PetscObject) z, "Reconstructed vector"));
120c4762a1bSJed Brown 
121c4762a1bSJed Brown       /* Set values of space vector x */
1225f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(x,rdm));
123c4762a1bSJed Brown 
1245f80ce2aSJacob Faibussowitsch       if (view) CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown       /* Apply FFTW_FORWARD and FFTW_BACKWARD */
1275f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(A,x,y));
1285f80ce2aSJacob Faibussowitsch       if (view) CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
129c4762a1bSJed Brown 
1305f80ce2aSJacob Faibussowitsch       CHKERRQ(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;
1345f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScale(z,a));
1355f80ce2aSJacob Faibussowitsch       if (view) CHKERRQ(VecView(z,PETSC_VIEWER_STDOUT_WORLD));
1365f80ce2aSJacob Faibussowitsch       CHKERRQ(VecAXPY(z,-1.0,x));
1375f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(z,NORM_1,&enorm));
138dd400576SPatrick Sanan       if (enorm > 1.e-9 && rank == 0) {
1395f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %e\n",enorm));
140c4762a1bSJed Brown       }
141c4762a1bSJed Brown 
1425f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&x));
1435f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&y));
1445f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&z));
1455f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&A));
146c4762a1bSJed Brown 
1475f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(dim));
148c4762a1bSJed Brown     }
149c4762a1bSJed Brown   }
150c4762a1bSJed Brown 
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rdm));
152*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
153*b122ec5aSJacob 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