xref: /petsc/src/mat/tests/ex121.c (revision ad0e0ea304bf4da9652c37aefec3cf21b1ef8ffc)
1c4762a1bSJed Brown static char help[] = "Test sequential FFTW convolution\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown   Compiling the code:
5c4762a1bSJed Brown     This code uses the complex numbers, so configure must be given --with-scalar-type=complex to enable this
6c4762a1bSJed Brown */
7c4762a1bSJed Brown 
8c4762a1bSJed Brown #include <petscmat.h>
9c4762a1bSJed Brown 
10*ad0e0ea3SStefano Zampini int main(int argc,char **args)
11c4762a1bSJed Brown {
12c4762a1bSJed Brown   typedef enum {RANDOM, CONSTANT, TANH, NUM_FUNCS} FuncType;
13c4762a1bSJed Brown   const char     *funcNames[NUM_FUNCS] = {"random", "constant", "tanh"};
14c4762a1bSJed Brown   Mat            A;
15c4762a1bSJed Brown   PetscMPIInt    size;
16c4762a1bSJed Brown   PetscInt       n = 10,N,ndim=4,dim[4],DIM,i,j;
17c4762a1bSJed Brown   Vec            w,x,y1,y2,z1,z2;
18c4762a1bSJed Brown   PetscScalar    *a, *a2, *a3;
19c4762a1bSJed Brown   PetscScalar    s;
20c4762a1bSJed Brown   PetscRandom    rdm;
21c4762a1bSJed Brown   PetscReal      enorm;
22c4762a1bSJed Brown   PetscInt       func     = 0;
23c4762a1bSJed Brown   FuncType       function = RANDOM;
24c4762a1bSJed Brown   PetscBool      view     = PETSC_FALSE;
25c4762a1bSJed Brown   PetscErrorCode ierr;
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
28c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRQ(ierr);
29c4762a1bSJed Brown   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This is a uniprocessor example only!");
30c4762a1bSJed Brown   ierr     = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex112");CHKERRQ(ierr);
31c4762a1bSJed Brown   ierr     = PetscOptionsEList("-function", "Function type", "ex121", funcNames, NUM_FUNCS, funcNames[function], &func, NULL);CHKERRQ(ierr);
32c4762a1bSJed Brown   ierr     = PetscOptionsBool("-vec_view draw", "View the functions", "ex112", view, &view, NULL);CHKERRQ(ierr);
33c4762a1bSJed Brown   function = (FuncType) func;
34c4762a1bSJed Brown   ierr     = PetscOptionsEnd();CHKERRQ(ierr);
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   for (DIM = 0; DIM < ndim; DIM++) {
37c4762a1bSJed Brown     dim[DIM] = n;  /* size of transformation in DIM-dimension */
38c4762a1bSJed Brown   }
39c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_SELF, &rdm);CHKERRQ(ierr);
40c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr);
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   for (DIM = 1; DIM < 5; DIM++) {
43c4762a1bSJed Brown     /* create vectors of length N=n^DIM */
44c4762a1bSJed Brown     for (i = 0, N = 1; i < DIM; i++) N *= dim[i];
45c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, "\n %d-D: FFTW on vector of size %d \n",DIM,N);CHKERRQ(ierr);
46c4762a1bSJed Brown     ierr = VecCreateSeq(PETSC_COMM_SELF,N,&x);CHKERRQ(ierr);
47c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) x, "Real space vector");CHKERRQ(ierr);
48c4762a1bSJed Brown     ierr = VecDuplicate(x,&w);CHKERRQ(ierr);
49c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) w, "Window vector");CHKERRQ(ierr);
50c4762a1bSJed Brown     ierr = VecDuplicate(x,&y1);CHKERRQ(ierr);
51c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) y1, "Frequency space vector");CHKERRQ(ierr);
52c4762a1bSJed Brown     ierr = VecDuplicate(x,&y2);CHKERRQ(ierr);
53c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) y2, "Frequency space window vector");CHKERRQ(ierr);
54c4762a1bSJed Brown     ierr = VecDuplicate(x,&z1);CHKERRQ(ierr);
55c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) z1, "Reconstructed convolution");CHKERRQ(ierr);
56c4762a1bSJed Brown     ierr = VecDuplicate(x,&z2);CHKERRQ(ierr);
57c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) z2, "Real space convolution");CHKERRQ(ierr);
58c4762a1bSJed Brown 
59c4762a1bSJed Brown     if (function == RANDOM) {
60c4762a1bSJed Brown       ierr = VecSetRandom(x, rdm);CHKERRQ(ierr);
61c4762a1bSJed Brown     } else if (function == CONSTANT) {
62c4762a1bSJed Brown       ierr = VecSet(x, 1.0);CHKERRQ(ierr);
63c4762a1bSJed Brown     } else if (function == TANH) {
64c4762a1bSJed Brown       ierr = VecGetArray(x, &a);CHKERRQ(ierr);
65c4762a1bSJed Brown       for (i = 0; i < N; ++i) {
66c4762a1bSJed Brown         a[i] = tanh((i - N/2.0)*(10.0/N));
67c4762a1bSJed Brown       }
68c4762a1bSJed Brown       ierr = VecRestoreArray(x, &a);CHKERRQ(ierr);
69c4762a1bSJed Brown     }
70c4762a1bSJed Brown     if (view) {ierr = VecView(x, PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);}
71c4762a1bSJed Brown 
72c4762a1bSJed Brown     /* Create window function */
73c4762a1bSJed Brown     ierr = VecGetArray(w, &a);CHKERRQ(ierr);
74c4762a1bSJed Brown     for (i = 0; i < N; ++i) {
75c4762a1bSJed Brown       /* Step Function */
76c4762a1bSJed Brown       a[i] = (i > N/4 && i < 3*N/4) ? 1.0 : 0.0;
77c4762a1bSJed Brown       /* Delta Function */
78c4762a1bSJed Brown       /*a[i] = (i == N/2)? 1.0: 0.0; */
79c4762a1bSJed Brown     }
80c4762a1bSJed Brown     ierr = VecRestoreArray(w, &a);CHKERRQ(ierr);
81c4762a1bSJed Brown     if (view) {ierr = VecView(w, PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);}
82c4762a1bSJed Brown 
83c4762a1bSJed Brown     /* create FFTW object */
84c4762a1bSJed Brown     ierr = MatCreateFFT(PETSC_COMM_SELF,DIM,dim,MATFFTW,&A);CHKERRQ(ierr);
85c4762a1bSJed Brown 
86c4762a1bSJed Brown     /* Convolve x with w*/
87c4762a1bSJed Brown     ierr = MatMult(A,x,y1);CHKERRQ(ierr);
88c4762a1bSJed Brown     ierr = MatMult(A,w,y2);CHKERRQ(ierr);
89c4762a1bSJed Brown     ierr = VecPointwiseMult(y1, y1, y2);CHKERRQ(ierr);
90c4762a1bSJed Brown     if (view && i == 0) {ierr = VecView(y1, PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);}
91c4762a1bSJed Brown     ierr = MatMultTranspose(A,y1,z1);CHKERRQ(ierr);
92c4762a1bSJed Brown 
93c4762a1bSJed Brown     /* Compute the real space convolution */
94c4762a1bSJed Brown     ierr = VecGetArray(x, &a);CHKERRQ(ierr);
95c4762a1bSJed Brown     ierr = VecGetArray(w, &a2);CHKERRQ(ierr);
96c4762a1bSJed Brown     ierr = VecGetArray(z2, &a3);CHKERRQ(ierr);
97c4762a1bSJed Brown     for (i = 0; i < N; ++i) {
98c4762a1bSJed Brown       /* PetscInt checkInd = (i > N/2-1)? i-N/2: i+N/2;*/
99c4762a1bSJed Brown 
100c4762a1bSJed Brown       a3[i] = 0.0;
101c4762a1bSJed Brown       for (j = -N/2+1; j < N/2; ++j) {
102c4762a1bSJed Brown         PetscInt xpInd   = (j < 0) ? N+j : j;
103c4762a1bSJed Brown         PetscInt diffInd = (i-j < 0) ? N-(j-i) : (i-j > N-1) ? i-j-N : i-j;
104c4762a1bSJed Brown 
105c4762a1bSJed Brown         a3[i] += a[xpInd]*a2[diffInd];
106c4762a1bSJed Brown       }
107c4762a1bSJed Brown     }
108c4762a1bSJed Brown     ierr = VecRestoreArray(x, &a);CHKERRQ(ierr);
109c4762a1bSJed Brown     ierr = VecRestoreArray(w, &a2);CHKERRQ(ierr);
110c4762a1bSJed Brown     ierr = VecRestoreArray(z2, &a3);CHKERRQ(ierr);
111c4762a1bSJed Brown 
112c4762a1bSJed Brown     /* compare z1 and z2. FFTW computes an unnormalized DFT, thus z1 = N*z2 */
113c4762a1bSJed Brown     s    = 1.0/(PetscReal)N;
114c4762a1bSJed Brown     ierr = VecScale(z1,s);CHKERRQ(ierr);
115c4762a1bSJed Brown     if (view) {ierr = VecView(z1, PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);}
116c4762a1bSJed Brown     if (view) {ierr = VecView(z2, PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);}
117c4762a1bSJed Brown     ierr = VecAXPY(z1,-1.0,z2);CHKERRQ(ierr);
118c4762a1bSJed Brown     ierr = VecNorm(z1,NORM_1,&enorm);CHKERRQ(ierr);
119c4762a1bSJed Brown     if (enorm > 1.e-11) {
120c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"  Error norm of |z1 - z2| %g\n",(double)enorm);CHKERRQ(ierr);
121c4762a1bSJed Brown     }
122c4762a1bSJed Brown 
123c4762a1bSJed Brown     /* free spaces */
124c4762a1bSJed Brown     ierr = VecDestroy(&x);CHKERRQ(ierr);
125c4762a1bSJed Brown     ierr = VecDestroy(&y1);CHKERRQ(ierr);
126c4762a1bSJed Brown     ierr = VecDestroy(&y2);CHKERRQ(ierr);
127c4762a1bSJed Brown     ierr = VecDestroy(&z1);CHKERRQ(ierr);
128c4762a1bSJed Brown     ierr = VecDestroy(&z2);CHKERRQ(ierr);
129c4762a1bSJed Brown     ierr = VecDestroy(&w);CHKERRQ(ierr);
130c4762a1bSJed Brown     ierr = MatDestroy(&A);CHKERRQ(ierr);
131c4762a1bSJed Brown   }
132c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr);
133c4762a1bSJed Brown   ierr = PetscFinalize();
134c4762a1bSJed Brown   return ierr;
135c4762a1bSJed Brown }
136c4762a1bSJed Brown 
137c4762a1bSJed Brown 
138c4762a1bSJed Brown /*TEST
139c4762a1bSJed Brown 
140c4762a1bSJed Brown    build:
141c4762a1bSJed Brown       requires: fftw complex
142c4762a1bSJed Brown 
143c4762a1bSJed Brown    test:
144c4762a1bSJed Brown       output_file: output/ex121.out
145c4762a1bSJed Brown       TODO: Example or FFTW interface is broken
146c4762a1bSJed Brown 
147c4762a1bSJed Brown TEST*/
148