xref: /petsc/src/mat/tests/ex121.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 
10ad0e0ea3SStefano 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 
27*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
285f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
292c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,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);
315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsEList("-function", "Function type", "ex121", funcNames, NUM_FUNCS, funcNames[function], &func, NULL));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-vec_view draw", "View the functions", "ex112", view, &view, NULL));
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   }
395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rdm));
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];
455f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "\n %d-D: FFTW on vector of size %d \n",DIM,N));
465f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,N,&x));
475f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) x, "Real space vector"));
485f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(x,&w));
495f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) w, "Window vector"));
505f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(x,&y1));
515f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) y1, "Frequency space vector"));
525f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(x,&y2));
535f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) y2, "Frequency space window vector"));
545f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(x,&z1));
555f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) z1, "Reconstructed convolution"));
565f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(x,&z2));
575f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) z2, "Real space convolution"));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown     if (function == RANDOM) {
605f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(x, rdm));
61c4762a1bSJed Brown     } else if (function == CONSTANT) {
625f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSet(x, 1.0));
63c4762a1bSJed Brown     } else if (function == TANH) {
645f80ce2aSJacob Faibussowitsch       CHKERRQ(VecGetArray(x, &a));
65c4762a1bSJed Brown       for (i = 0; i < N; ++i) {
66c4762a1bSJed Brown         a[i] = tanh((i - N/2.0)*(10.0/N));
67c4762a1bSJed Brown       }
685f80ce2aSJacob Faibussowitsch       CHKERRQ(VecRestoreArray(x, &a));
69c4762a1bSJed Brown     }
705f80ce2aSJacob Faibussowitsch     if (view) CHKERRQ(VecView(x, PETSC_VIEWER_DRAW_WORLD));
71c4762a1bSJed Brown 
72c4762a1bSJed Brown     /* Create window function */
735f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(w, &a));
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     }
805f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(w, &a));
815f80ce2aSJacob Faibussowitsch     if (view) CHKERRQ(VecView(w, PETSC_VIEWER_DRAW_WORLD));
82c4762a1bSJed Brown 
83c4762a1bSJed Brown     /* create FFTW object */
845f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateFFT(PETSC_COMM_SELF,DIM,dim,MATFFTW,&A));
85c4762a1bSJed Brown 
86c4762a1bSJed Brown     /* Convolve x with w*/
875f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(A,x,y1));
885f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(A,w,y2));
895f80ce2aSJacob Faibussowitsch     CHKERRQ(VecPointwiseMult(y1, y1, y2));
905f80ce2aSJacob Faibussowitsch     if (view && i == 0) CHKERRQ(VecView(y1, PETSC_VIEWER_DRAW_WORLD));
915f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultTranspose(A,y1,z1));
92c4762a1bSJed Brown 
93c4762a1bSJed Brown     /* Compute the real space convolution */
945f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(x, &a));
955f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(w, &a2));
965f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(z2, &a3));
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     }
1085f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(x, &a));
1095f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(w, &a2));
1105f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(z2, &a3));
111c4762a1bSJed Brown 
112c4762a1bSJed Brown     /* compare z1 and z2. FFTW computes an unnormalized DFT, thus z1 = N*z2 */
113c4762a1bSJed Brown     s    = 1.0/(PetscReal)N;
1145f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScale(z1,s));
1155f80ce2aSJacob Faibussowitsch     if (view) CHKERRQ(VecView(z1, PETSC_VIEWER_DRAW_WORLD));
1165f80ce2aSJacob Faibussowitsch     if (view) CHKERRQ(VecView(z2, PETSC_VIEWER_DRAW_WORLD));
1175f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(z1,-1.0,z2));
1185f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(z1,NORM_1,&enorm));
119c4762a1bSJed Brown     if (enorm > 1.e-11) {
1205f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"  Error norm of |z1 - z2| %g\n",(double)enorm));
121c4762a1bSJed Brown     }
122c4762a1bSJed Brown 
123c4762a1bSJed Brown     /* free spaces */
1245f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&x));
1255f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&y1));
1265f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&y2));
1275f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&z1));
1285f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&z2));
1295f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&w));
1305f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&A));
131c4762a1bSJed Brown   }
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rdm));
133*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
134*b122ec5aSJacob Faibussowitsch   return 0;
135c4762a1bSJed Brown }
136c4762a1bSJed Brown 
137c4762a1bSJed Brown /*TEST
138c4762a1bSJed Brown 
139c4762a1bSJed Brown    build:
140c4762a1bSJed Brown       requires: fftw complex
141c4762a1bSJed Brown 
142c4762a1bSJed Brown    test:
143c4762a1bSJed Brown       output_file: output/ex121.out
144c4762a1bSJed Brown       TODO: Example or FFTW interface is broken
145c4762a1bSJed Brown 
146c4762a1bSJed Brown TEST*/
147