xref: /petsc/src/mat/tests/ex143.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Illustrate how to use mpi FFTW and PETSc-FFTW interface \n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown /*
4*c4762a1bSJed Brown   Compiling the code:
5*c4762a1bSJed Brown       This code uses the complex numbers version of PETSc, so configure
6*c4762a1bSJed Brown       must be run to enable this
7*c4762a1bSJed Brown 
8*c4762a1bSJed Brown  Usage:
9*c4762a1bSJed Brown    mpiexec -n <np> ./ex143 -use_FFTW_interface NO
10*c4762a1bSJed Brown    mpiexec -n <np> ./ex143 -use_FFTW_interface YES
11*c4762a1bSJed Brown */
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown #include <petscmat.h>
14*c4762a1bSJed Brown #include <fftw3-mpi.h>
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown int main(int argc,char **args)
17*c4762a1bSJed Brown {
18*c4762a1bSJed Brown   PetscErrorCode ierr;
19*c4762a1bSJed Brown   PetscMPIInt    rank,size;
20*c4762a1bSJed Brown   PetscInt       N0=50,N1=20,N=N0*N1,DIM;
21*c4762a1bSJed Brown   PetscRandom    rdm;
22*c4762a1bSJed Brown   PetscScalar    a;
23*c4762a1bSJed Brown   PetscReal      enorm;
24*c4762a1bSJed Brown   Vec            x,y,z;
25*c4762a1bSJed Brown   PetscBool      view=PETSC_FALSE,use_interface=PETSC_TRUE;
26*c4762a1bSJed Brown 
27*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
28*c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex143");CHKERRQ(ierr);
29*c4762a1bSJed Brown   ierr = PetscOptionsBool("-vec_view draw", "View the vectors", "ex143", view, &view, NULL);CHKERRQ(ierr);
30*c4762a1bSJed Brown   ierr = PetscOptionsBool("-use_FFTW_interface", "Use PETSc-FFTW interface", "ex143",use_interface, &use_interface, NULL);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-use_FFTW_interface",&use_interface,NULL);CHKERRQ(ierr);
34*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRQ(ierr);
35*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr);
36*c4762a1bSJed Brown 
37*c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_WORLD, &rdm);CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr);
39*c4762a1bSJed Brown 
40*c4762a1bSJed Brown   if (!use_interface) {
41*c4762a1bSJed Brown     /* Use mpi FFTW without PETSc-FFTW interface, 2D case only */
42*c4762a1bSJed Brown     /*---------------------------------------------------------*/
43*c4762a1bSJed Brown     fftw_plan    fplan,bplan;
44*c4762a1bSJed Brown     fftw_complex *data_in,*data_out,*data_out2;
45*c4762a1bSJed Brown     ptrdiff_t    alloc_local,local_n0,local_0_start;
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown     DIM = 2;
48*c4762a1bSJed Brown     if (!rank) {
49*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"Use FFTW without PETSc-FFTW interface, DIM %D\n",DIM);CHKERRQ(ierr);
50*c4762a1bSJed Brown     }
51*c4762a1bSJed Brown     fftw_mpi_init();
52*c4762a1bSJed Brown     N           = N0*N1;
53*c4762a1bSJed Brown     alloc_local = fftw_mpi_local_size_2d(N0,N1,PETSC_COMM_WORLD,&local_n0,&local_0_start);
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown     data_in   = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
56*c4762a1bSJed Brown     data_out  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
57*c4762a1bSJed Brown     data_out2 = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
58*c4762a1bSJed Brown 
59*c4762a1bSJed Brown     ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_in,&x);CHKERRQ(ierr);
60*c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) x, "Real Space vector");CHKERRQ(ierr);
61*c4762a1bSJed Brown     ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out,&y);CHKERRQ(ierr);
62*c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) y, "Frequency space vector");CHKERRQ(ierr);
63*c4762a1bSJed Brown     ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,(PetscInt)local_n0*N1,(PetscInt)N,(const PetscScalar*)data_out2,&z);CHKERRQ(ierr);
64*c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) z, "Reconstructed vector");CHKERRQ(ierr);
65*c4762a1bSJed Brown 
66*c4762a1bSJed Brown     fplan = fftw_mpi_plan_dft_2d(N0,N1,data_in,data_out,PETSC_COMM_WORLD,FFTW_FORWARD,FFTW_ESTIMATE);
67*c4762a1bSJed Brown     bplan = fftw_mpi_plan_dft_2d(N0,N1,data_out,data_out2,PETSC_COMM_WORLD,FFTW_BACKWARD,FFTW_ESTIMATE);
68*c4762a1bSJed Brown 
69*c4762a1bSJed Brown     ierr = VecSetRandom(x, rdm);CHKERRQ(ierr);
70*c4762a1bSJed Brown     if (view) {ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);}
71*c4762a1bSJed Brown 
72*c4762a1bSJed Brown     fftw_execute(fplan);
73*c4762a1bSJed Brown     if (view) {ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);}
74*c4762a1bSJed Brown 
75*c4762a1bSJed Brown     fftw_execute(bplan);
76*c4762a1bSJed Brown 
77*c4762a1bSJed Brown     /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */
78*c4762a1bSJed Brown     a    = 1.0/(PetscReal)N;
79*c4762a1bSJed Brown     ierr = VecScale(z,a);CHKERRQ(ierr);
80*c4762a1bSJed Brown     if (view) {ierr = VecView(z, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);}
81*c4762a1bSJed Brown     ierr = VecAXPY(z,-1.0,x);CHKERRQ(ierr);
82*c4762a1bSJed Brown     ierr = VecNorm(z,NORM_1,&enorm);CHKERRQ(ierr);
83*c4762a1bSJed Brown     if (enorm > 1.e-11 && !rank) {
84*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %g\n",(double)enorm);CHKERRQ(ierr);
85*c4762a1bSJed Brown     }
86*c4762a1bSJed Brown 
87*c4762a1bSJed Brown     /* Free spaces */
88*c4762a1bSJed Brown     fftw_destroy_plan(fplan);
89*c4762a1bSJed Brown     fftw_destroy_plan(bplan);
90*c4762a1bSJed Brown     fftw_free(data_in);  ierr = VecDestroy(&x);CHKERRQ(ierr);
91*c4762a1bSJed Brown     fftw_free(data_out); ierr = VecDestroy(&y);CHKERRQ(ierr);
92*c4762a1bSJed Brown     fftw_free(data_out2);ierr = VecDestroy(&z);CHKERRQ(ierr);
93*c4762a1bSJed Brown 
94*c4762a1bSJed Brown   } else {
95*c4762a1bSJed Brown     /* Use PETSc-FFTW interface                  */
96*c4762a1bSJed Brown     /*-------------------------------------------*/
97*c4762a1bSJed Brown     PetscInt i,*dim,k;
98*c4762a1bSJed Brown     Mat      A;
99*c4762a1bSJed Brown 
100*c4762a1bSJed Brown     N=1;
101*c4762a1bSJed Brown     for (i=1; i<5; i++) {
102*c4762a1bSJed Brown       DIM  = i;
103*c4762a1bSJed Brown       ierr = PetscMalloc1(i,&dim);CHKERRQ(ierr);
104*c4762a1bSJed Brown       for (k=0; k<i; k++) {
105*c4762a1bSJed Brown         dim[k]=30;
106*c4762a1bSJed Brown       }
107*c4762a1bSJed Brown       N *= dim[i-1];
108*c4762a1bSJed Brown 
109*c4762a1bSJed Brown 
110*c4762a1bSJed Brown       /* Create FFTW object */
111*c4762a1bSJed Brown       if (!rank) printf("Use PETSc-FFTW interface...%d-DIM: %d\n",(int)DIM,(int)N);
112*c4762a1bSJed Brown 
113*c4762a1bSJed Brown       ierr = MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A);CHKERRQ(ierr);
114*c4762a1bSJed Brown 
115*c4762a1bSJed Brown       /* Create vectors that are compatible with parallel layout of A - must call MatCreateVecs()! */
116*c4762a1bSJed Brown 
117*c4762a1bSJed Brown       ierr = MatCreateVecsFFTW(A,&x,&y,&z);CHKERRQ(ierr);
118*c4762a1bSJed Brown       ierr = PetscObjectSetName((PetscObject) x, "Real space vector");CHKERRQ(ierr);
119*c4762a1bSJed Brown       ierr = PetscObjectSetName((PetscObject) y, "Frequency space vector");CHKERRQ(ierr);
120*c4762a1bSJed Brown       ierr = PetscObjectSetName((PetscObject) z, "Reconstructed vector");CHKERRQ(ierr);
121*c4762a1bSJed Brown 
122*c4762a1bSJed Brown       /* Set values of space vector x */
123*c4762a1bSJed Brown       ierr = VecSetRandom(x,rdm);CHKERRQ(ierr);
124*c4762a1bSJed Brown 
125*c4762a1bSJed Brown       if (view) {ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);}
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown       /* Apply FFTW_FORWARD and FFTW_BACKWARD */
128*c4762a1bSJed Brown       ierr = MatMult(A,x,y);CHKERRQ(ierr);
129*c4762a1bSJed Brown       if (view) {ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);}
130*c4762a1bSJed Brown 
131*c4762a1bSJed Brown       ierr = MatMultTranspose(A,y,z);CHKERRQ(ierr);
132*c4762a1bSJed Brown 
133*c4762a1bSJed Brown       /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */
134*c4762a1bSJed Brown       a    = 1.0/(PetscReal)N;
135*c4762a1bSJed Brown       ierr = VecScale(z,a);CHKERRQ(ierr);
136*c4762a1bSJed Brown       if (view) {ierr = VecView(z,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);}
137*c4762a1bSJed Brown       ierr = VecAXPY(z,-1.0,x);CHKERRQ(ierr);
138*c4762a1bSJed Brown       ierr = VecNorm(z,NORM_1,&enorm);CHKERRQ(ierr);
139*c4762a1bSJed Brown       if (enorm > 1.e-9 && !rank) {
140*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %e\n",enorm);CHKERRQ(ierr);
141*c4762a1bSJed Brown       }
142*c4762a1bSJed Brown 
143*c4762a1bSJed Brown       ierr = VecDestroy(&x);CHKERRQ(ierr);
144*c4762a1bSJed Brown       ierr = VecDestroy(&y);CHKERRQ(ierr);
145*c4762a1bSJed Brown       ierr = VecDestroy(&z);CHKERRQ(ierr);
146*c4762a1bSJed Brown       ierr = MatDestroy(&A);CHKERRQ(ierr);
147*c4762a1bSJed Brown 
148*c4762a1bSJed Brown       ierr = PetscFree(dim);CHKERRQ(ierr);
149*c4762a1bSJed Brown     }
150*c4762a1bSJed Brown   }
151*c4762a1bSJed Brown 
152*c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr);
153*c4762a1bSJed Brown   ierr = PetscFinalize();
154*c4762a1bSJed Brown   return ierr;
155*c4762a1bSJed Brown }
156*c4762a1bSJed Brown 
157*c4762a1bSJed Brown 
158*c4762a1bSJed Brown /*TEST
159*c4762a1bSJed Brown 
160*c4762a1bSJed Brown    build:
161*c4762a1bSJed Brown       requires: fftw complex
162*c4762a1bSJed Brown 
163*c4762a1bSJed Brown    test:
164*c4762a1bSJed Brown       output_file: output/ex143.out
165*c4762a1bSJed Brown 
166*c4762a1bSJed Brown    test:
167*c4762a1bSJed Brown       suffix: 2
168*c4762a1bSJed Brown       nsize: 3
169*c4762a1bSJed Brown       output_file: output/ex143.out
170*c4762a1bSJed Brown 
171*c4762a1bSJed Brown TEST*/
172