xref: /petsc/src/mat/tests/ex142.c (revision 030f984af8d8bb4c203755d35bded3c05b3d83ce)
1 static char help[] = "Test sequential r2c/c2r FFTW without PETSc interface \n\n";
2 
3 /*
4   Compiling the code:
5       This code uses the real numbers version of PETSc
6 */
7 
8 #include <petscmat.h>
9 #include <fftw3.h>
10 
11 int main(int argc,char **args)
12 {
13   typedef enum {RANDOM, CONSTANT, TANH, NUM_FUNCS} FuncType;
14   const char      *funcNames[NUM_FUNCS] = {"random", "constant", "tanh"};
15   PetscMPIInt     size;
16   int             n = 10,N,Ny,ndim=4,i,dim[4],DIM;
17   Vec             x,y,z;
18   PetscScalar     s;
19   PetscRandom     rdm;
20   PetscReal       enorm;
21   PetscInt        func     = RANDOM;
22   FuncType        function = RANDOM;
23   PetscBool       view     = PETSC_FALSE;
24   PetscErrorCode  ierr;
25   PetscScalar     *x_array,*y_array,*z_array;
26   fftw_plan       fplan,bplan;
27 
28   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
29 #if defined(PETSC_USE_COMPLEX)
30   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers");
31 #endif
32 
33   ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRMPI(ierr);
34   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This is a uniprocessor example only!");
35   ierr     = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex142");CHKERRQ(ierr);
36   ierr     = PetscOptionsEList("-function", "Function type", "ex142", funcNames, NUM_FUNCS, funcNames[function], &func, NULL);CHKERRQ(ierr);
37   ierr     = PetscOptionsBool("-vec_view draw", "View the functions", "ex142", view, &view, NULL);CHKERRQ(ierr);
38   function = (FuncType) func;
39   ierr     = PetscOptionsEnd();CHKERRQ(ierr);
40 
41   for (DIM = 0; DIM < ndim; DIM++) {
42     dim[DIM] = n;  /* size of real space vector in DIM-dimension */
43   }
44   ierr = PetscRandomCreate(PETSC_COMM_SELF, &rdm);CHKERRQ(ierr);
45   ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr);
46 
47   for (DIM = 1; DIM < 5; DIM++) {
48     /* create vectors of length N=dim[0]*dim[1]* ...*dim[DIM-1] */
49     /*----------------------------------------------------------*/
50     N = Ny = 1;
51     for (i = 0; i < DIM-1; i++) {
52       N *= dim[i];
53     }
54     Ny = N; Ny *= 2*(dim[DIM-1]/2 + 1); /* add padding elements to output vector y */
55     N *= dim[DIM-1];
56 
57     ierr = PetscPrintf(PETSC_COMM_SELF, "\n %d-D: FFTW on vector of size %d \n",DIM,N);CHKERRQ(ierr);
58     ierr = VecCreateSeq(PETSC_COMM_SELF,N,&x);CHKERRQ(ierr);
59     ierr = PetscObjectSetName((PetscObject) x, "Real space vector");CHKERRQ(ierr);
60 
61     ierr = VecCreateSeq(PETSC_COMM_SELF,Ny,&y);CHKERRQ(ierr);
62     ierr = PetscObjectSetName((PetscObject) y, "Frequency space vector");CHKERRQ(ierr);
63 
64     ierr = VecDuplicate(x,&z);CHKERRQ(ierr);
65     ierr = PetscObjectSetName((PetscObject) z, "Reconstructed vector");CHKERRQ(ierr);
66 
67     /* Set fftw plan                    */
68     /*----------------------------------*/
69     ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
70     ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
71     ierr = VecGetArray(z,&z_array);CHKERRQ(ierr);
72 
73     unsigned int flags = FFTW_ESTIMATE; /*or FFTW_MEASURE */
74     /* The data in the in/out arrays is overwritten during FFTW_MEASURE planning, so such planning
75      should be done before the input is initialized by the user. */
76     ierr = PetscPrintf(PETSC_COMM_SELF,"DIM: %d, N %d, Ny %d\n",DIM,N,Ny);CHKERRQ(ierr);
77 
78     switch (DIM) {
79     case 1:
80       fplan = fftw_plan_dft_r2c_1d(dim[0], (double*)x_array, (fftw_complex*)y_array, flags);
81       bplan = fftw_plan_dft_c2r_1d(dim[0], (fftw_complex*)y_array, (double*)z_array, flags);
82       break;
83     case 2:
84       fplan = fftw_plan_dft_r2c_2d(dim[0],dim[1],(double*)x_array, (fftw_complex*)y_array,flags);
85       bplan = fftw_plan_dft_c2r_2d(dim[0],dim[1],(fftw_complex*)y_array,(double*)z_array,flags);
86       break;
87     case 3:
88       fplan = fftw_plan_dft_r2c_3d(dim[0],dim[1],dim[2],(double*)x_array, (fftw_complex*)y_array,flags);
89       bplan = fftw_plan_dft_c2r_3d(dim[0],dim[1],dim[2],(fftw_complex*)y_array,(double*)z_array,flags);
90       break;
91     default:
92       fplan = fftw_plan_dft_r2c(DIM,(int*)dim,(double*)x_array, (fftw_complex*)y_array,flags);
93       bplan = fftw_plan_dft_c2r(DIM,(int*)dim,(fftw_complex*)y_array,(double*)z_array,flags);
94       break;
95     }
96 
97     ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
98     ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
99     ierr = VecRestoreArray(z,&z_array);CHKERRQ(ierr);
100 
101     /* Initialize Real space vector x:
102        The data in the in/out arrays is overwritten during FFTW_MEASURE planning, so planning
103        should be done before the input is initialized by the user.
104     --------------------------------------------------------*/
105     if (function == RANDOM) {
106       ierr = VecSetRandom(x, rdm);CHKERRQ(ierr);
107     } else if (function == CONSTANT) {
108       ierr = VecSet(x, 1.0);CHKERRQ(ierr);
109     } else if (function == TANH) {
110       ierr = VecGetArray(x, &x_array);CHKERRQ(ierr);
111       for (i = 0; i < N; ++i) {
112         x_array[i] = tanh((i - N/2.0)*(10.0/N));
113       }
114       ierr = VecRestoreArray(x, &x_array);CHKERRQ(ierr);
115     }
116     if (view) {
117       ierr = VecView(x, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
118     }
119 
120     /* FFT - also test repeated transformation   */
121     /*-------------------------------------------*/
122     ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
123     ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
124     ierr = VecGetArray(z,&z_array);CHKERRQ(ierr);
125     for (i=0; i<4; i++) {
126       /* FFTW_FORWARD */
127       fftw_execute(fplan);
128 
129       /* FFTW_BACKWARD: destroys its input array 'y_array' even for out-of-place transforms! */
130       fftw_execute(bplan);
131     }
132     ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
133     ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
134     ierr = VecRestoreArray(z,&z_array);CHKERRQ(ierr);
135 
136     /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */
137     /*------------------------------------------------------------------*/
138     s    = 1.0/(PetscReal)N;
139     ierr = VecScale(z,s);CHKERRQ(ierr);
140     if (view) {ierr = VecView(x, PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);}
141     if (view) {ierr = VecView(z, PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);}
142     ierr = VecAXPY(z,-1.0,x);CHKERRQ(ierr);
143     ierr = VecNorm(z,NORM_1,&enorm);CHKERRQ(ierr);
144     if (enorm > 1.e-11) {
145       ierr = PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %g\n",(double)enorm);CHKERRQ(ierr);
146     }
147 
148     /* free spaces */
149     fftw_destroy_plan(fplan);
150     fftw_destroy_plan(bplan);
151     ierr = VecDestroy(&x);CHKERRQ(ierr);
152     ierr = VecDestroy(&y);CHKERRQ(ierr);
153     ierr = VecDestroy(&z);CHKERRQ(ierr);
154   }
155   ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr);
156   ierr = PetscFinalize();
157   return ierr;
158 }
159 
160 /*TEST
161 
162    build:
163      requires: fftw !complex
164 
165    test:
166      output_file: output/ex142.out
167 
168 TEST*/
169