xref: /petsc/src/mat/tests/ex121.c (revision d0609ced746bc51b019815ca91d747429db24893)
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 
269566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
28be096a46SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
29*d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex112");
309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-function", "Function type", "ex121", funcNames, NUM_FUNCS, funcNames[function], &func, NULL));
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-vec_view draw", "View the functions", "ex112", view, &view, NULL));
32c4762a1bSJed Brown   function = (FuncType) func;
33*d0609cedSBarry Smith   PetscOptionsEnd();
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   for (DIM = 0; DIM < ndim; DIM++) {
36c4762a1bSJed Brown     dim[DIM] = n;  /* size of transformation in DIM-dimension */
37c4762a1bSJed Brown   }
389566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
399566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   for (DIM = 1; DIM < 5; DIM++) {
42c4762a1bSJed Brown     /* create vectors of length N=n^DIM */
43c4762a1bSJed Brown     for (i = 0, N = 1; i < DIM; i++) N *= dim[i];
449566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n %d-D: FFTW on vector of size %d \n",DIM,N));
459566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF,N,&x));
469566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) x, "Real space vector"));
479566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x,&w));
489566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) w, "Window vector"));
499566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x,&y1));
509566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) y1, "Frequency space vector"));
519566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x,&y2));
529566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) y2, "Frequency space window vector"));
539566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x,&z1));
549566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) z1, "Reconstructed convolution"));
559566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x,&z2));
569566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) z2, "Real space convolution"));
57c4762a1bSJed Brown 
58c4762a1bSJed Brown     if (function == RANDOM) {
599566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(x, rdm));
60c4762a1bSJed Brown     } else if (function == CONSTANT) {
619566063dSJacob Faibussowitsch       PetscCall(VecSet(x, 1.0));
62c4762a1bSJed Brown     } else if (function == TANH) {
639566063dSJacob Faibussowitsch       PetscCall(VecGetArray(x, &a));
64c4762a1bSJed Brown       for (i = 0; i < N; ++i) {
65c4762a1bSJed Brown         a[i] = tanh((i - N/2.0)*(10.0/N));
66c4762a1bSJed Brown       }
679566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(x, &a));
68c4762a1bSJed Brown     }
699566063dSJacob Faibussowitsch     if (view) PetscCall(VecView(x, PETSC_VIEWER_DRAW_WORLD));
70c4762a1bSJed Brown 
71c4762a1bSJed Brown     /* Create window function */
729566063dSJacob Faibussowitsch     PetscCall(VecGetArray(w, &a));
73c4762a1bSJed Brown     for (i = 0; i < N; ++i) {
74c4762a1bSJed Brown       /* Step Function */
75c4762a1bSJed Brown       a[i] = (i > N/4 && i < 3*N/4) ? 1.0 : 0.0;
76c4762a1bSJed Brown       /* Delta Function */
77c4762a1bSJed Brown       /*a[i] = (i == N/2)? 1.0: 0.0; */
78c4762a1bSJed Brown     }
799566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(w, &a));
809566063dSJacob Faibussowitsch     if (view) PetscCall(VecView(w, PETSC_VIEWER_DRAW_WORLD));
81c4762a1bSJed Brown 
82c4762a1bSJed Brown     /* create FFTW object */
839566063dSJacob Faibussowitsch     PetscCall(MatCreateFFT(PETSC_COMM_SELF,DIM,dim,MATFFTW,&A));
84c4762a1bSJed Brown 
85c4762a1bSJed Brown     /* Convolve x with w*/
869566063dSJacob Faibussowitsch     PetscCall(MatMult(A,x,y1));
879566063dSJacob Faibussowitsch     PetscCall(MatMult(A,w,y2));
889566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(y1, y1, y2));
899566063dSJacob Faibussowitsch     if (view && i == 0) PetscCall(VecView(y1, PETSC_VIEWER_DRAW_WORLD));
909566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A,y1,z1));
91c4762a1bSJed Brown 
92c4762a1bSJed Brown     /* Compute the real space convolution */
939566063dSJacob Faibussowitsch     PetscCall(VecGetArray(x, &a));
949566063dSJacob Faibussowitsch     PetscCall(VecGetArray(w, &a2));
959566063dSJacob Faibussowitsch     PetscCall(VecGetArray(z2, &a3));
96c4762a1bSJed Brown     for (i = 0; i < N; ++i) {
97c4762a1bSJed Brown       /* PetscInt checkInd = (i > N/2-1)? i-N/2: i+N/2;*/
98c4762a1bSJed Brown 
99c4762a1bSJed Brown       a3[i] = 0.0;
100c4762a1bSJed Brown       for (j = -N/2+1; j < N/2; ++j) {
101c4762a1bSJed Brown         PetscInt xpInd   = (j < 0) ? N+j : j;
102c4762a1bSJed Brown         PetscInt diffInd = (i-j < 0) ? N-(j-i) : (i-j > N-1) ? i-j-N : i-j;
103c4762a1bSJed Brown 
104c4762a1bSJed Brown         a3[i] += a[xpInd]*a2[diffInd];
105c4762a1bSJed Brown       }
106c4762a1bSJed Brown     }
1079566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(x, &a));
1089566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(w, &a2));
1099566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(z2, &a3));
110c4762a1bSJed Brown 
111c4762a1bSJed Brown     /* compare z1 and z2. FFTW computes an unnormalized DFT, thus z1 = N*z2 */
112c4762a1bSJed Brown     s    = 1.0/(PetscReal)N;
1139566063dSJacob Faibussowitsch     PetscCall(VecScale(z1,s));
1149566063dSJacob Faibussowitsch     if (view) PetscCall(VecView(z1, PETSC_VIEWER_DRAW_WORLD));
1159566063dSJacob Faibussowitsch     if (view) PetscCall(VecView(z2, PETSC_VIEWER_DRAW_WORLD));
1169566063dSJacob Faibussowitsch     PetscCall(VecAXPY(z1,-1.0,z2));
1179566063dSJacob Faibussowitsch     PetscCall(VecNorm(z1,NORM_1,&enorm));
118c4762a1bSJed Brown     if (enorm > 1.e-11) {
1199566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF,"  Error norm of |z1 - z2| %g\n",(double)enorm));
120c4762a1bSJed Brown     }
121c4762a1bSJed Brown 
122c4762a1bSJed Brown     /* free spaces */
1239566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&x));
1249566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&y1));
1259566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&y2));
1269566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&z1));
1279566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&z2));
1289566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&w));
1299566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
130c4762a1bSJed Brown   }
1319566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
1329566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
133b122ec5aSJacob Faibussowitsch   return 0;
134c4762a1bSJed Brown }
135c4762a1bSJed Brown 
136c4762a1bSJed Brown /*TEST
137c4762a1bSJed Brown 
138c4762a1bSJed Brown    build:
139c4762a1bSJed Brown       requires: fftw complex
140c4762a1bSJed Brown 
141c4762a1bSJed Brown    test:
142c4762a1bSJed Brown       output_file: output/ex121.out
143c4762a1bSJed Brown       TODO: Example or FFTW interface is broken
144c4762a1bSJed Brown 
145c4762a1bSJed Brown TEST*/
146