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