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
main(int argc,char ** args)10d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
11d71ae5a4SJacob Faibussowitsch {
129371c9d4SSatish Balay typedef enum {
139371c9d4SSatish Balay RANDOM,
149371c9d4SSatish Balay CONSTANT,
159371c9d4SSatish Balay TANH,
169371c9d4SSatish Balay NUM_FUNCS
179371c9d4SSatish Balay } FuncType;
18c4762a1bSJed Brown const char *funcNames[NUM_FUNCS] = {"random", "constant", "tanh"};
19c4762a1bSJed Brown Mat A;
20c4762a1bSJed Brown PetscMPIInt size;
21c4762a1bSJed Brown PetscInt n = 10, N, ndim = 4, dim[4], DIM, i, j;
22c4762a1bSJed Brown Vec w, x, y1, y2, z1, z2;
23c4762a1bSJed Brown PetscScalar *a, *a2, *a3;
24c4762a1bSJed Brown PetscScalar s;
25c4762a1bSJed Brown PetscRandom rdm;
26c4762a1bSJed Brown PetscReal enorm;
27c4762a1bSJed Brown PetscInt func = 0;
28c4762a1bSJed Brown FuncType function = RANDOM;
29c4762a1bSJed Brown PetscBool view = PETSC_FALSE;
30c4762a1bSJed Brown
31327415f7SBarry Smith PetscFunctionBeginUser;
32c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
34be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
35d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex112");
369566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-function", "Function type", "ex121", funcNames, NUM_FUNCS, funcNames[function], &func, NULL));
379566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-vec_view draw", "View the functions", "ex112", view, &view, NULL));
38c4762a1bSJed Brown function = (FuncType)func;
39d0609cedSBarry Smith PetscOptionsEnd();
40c4762a1bSJed Brown
41*ac530a7eSPierre Jolivet for (DIM = 0; DIM < ndim; DIM++) dim[DIM] = n; /* size of transformation in DIM-dimension */
429566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
439566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm));
44c4762a1bSJed Brown
45c4762a1bSJed Brown for (DIM = 1; DIM < 5; DIM++) {
46c4762a1bSJed Brown /* create vectors of length N=n^DIM */
47c4762a1bSJed Brown for (i = 0, N = 1; i < DIM; i++) N *= dim[i];
489566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n %d-D: FFTW on vector of size %d \n", DIM, N));
499566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, N, &x));
509566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x, "Real space vector"));
519566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &w));
529566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)w, "Window vector"));
539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y1));
549566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)y1, "Frequency space vector"));
559566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y2));
569566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)y2, "Frequency space window vector"));
579566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &z1));
589566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)z1, "Reconstructed convolution"));
599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &z2));
609566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)z2, "Real space convolution"));
61c4762a1bSJed Brown
62c4762a1bSJed Brown if (function == RANDOM) {
639566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm));
64c4762a1bSJed Brown } else if (function == CONSTANT) {
659566063dSJacob Faibussowitsch PetscCall(VecSet(x, 1.0));
66c4762a1bSJed Brown } else if (function == TANH) {
679566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &a));
68ad540459SPierre Jolivet for (i = 0; i < N; ++i) a[i] = tanh((i - N / 2.0) * (10.0 / N));
699566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &a));
70c4762a1bSJed Brown }
719566063dSJacob Faibussowitsch if (view) PetscCall(VecView(x, PETSC_VIEWER_DRAW_WORLD));
72c4762a1bSJed Brown
73c4762a1bSJed Brown /* Create window function */
749566063dSJacob Faibussowitsch PetscCall(VecGetArray(w, &a));
75c4762a1bSJed Brown for (i = 0; i < N; ++i) {
76c4762a1bSJed Brown /* Step Function */
77c4762a1bSJed Brown a[i] = (i > N / 4 && i < 3 * N / 4) ? 1.0 : 0.0;
78c4762a1bSJed Brown /* Delta Function */
79c4762a1bSJed Brown /*a[i] = (i == N/2)? 1.0: 0.0; */
80c4762a1bSJed Brown }
819566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w, &a));
829566063dSJacob Faibussowitsch if (view) PetscCall(VecView(w, PETSC_VIEWER_DRAW_WORLD));
83c4762a1bSJed Brown
84c4762a1bSJed Brown /* create FFTW object */
859566063dSJacob Faibussowitsch PetscCall(MatCreateFFT(PETSC_COMM_SELF, DIM, dim, MATFFTW, &A));
86c4762a1bSJed Brown
87c4762a1bSJed Brown /* Convolve x with w*/
889566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, y1));
899566063dSJacob Faibussowitsch PetscCall(MatMult(A, w, y2));
909566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(y1, y1, y2));
919566063dSJacob Faibussowitsch if (view && i == 0) PetscCall(VecView(y1, PETSC_VIEWER_DRAW_WORLD));
929566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, y1, z1));
93c4762a1bSJed Brown
94c4762a1bSJed Brown /* Compute the real space convolution */
959566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &a));
969566063dSJacob Faibussowitsch PetscCall(VecGetArray(w, &a2));
979566063dSJacob Faibussowitsch PetscCall(VecGetArray(z2, &a3));
98c4762a1bSJed Brown for (i = 0; i < N; ++i) {
99c4762a1bSJed Brown /* PetscInt checkInd = (i > N/2-1)? i-N/2: i+N/2;*/
100c4762a1bSJed Brown
101c4762a1bSJed Brown a3[i] = 0.0;
102c4762a1bSJed Brown for (j = -N / 2 + 1; j < N / 2; ++j) {
103c4762a1bSJed Brown PetscInt xpInd = (j < 0) ? N + j : j;
104c4762a1bSJed Brown PetscInt diffInd = (i - j < 0) ? N - (j - i) : (i - j > N - 1) ? i - j - N : i - j;
105c4762a1bSJed Brown
106c4762a1bSJed Brown a3[i] += a[xpInd] * a2[diffInd];
107c4762a1bSJed Brown }
108c4762a1bSJed Brown }
1099566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &a));
1109566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w, &a2));
1119566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(z2, &a3));
112c4762a1bSJed Brown
113c4762a1bSJed Brown /* compare z1 and z2. FFTW computes an unnormalized DFT, thus z1 = N*z2 */
114c4762a1bSJed Brown s = 1.0 / (PetscReal)N;
1159566063dSJacob Faibussowitsch PetscCall(VecScale(z1, s));
1169566063dSJacob Faibussowitsch if (view) PetscCall(VecView(z1, PETSC_VIEWER_DRAW_WORLD));
1179566063dSJacob Faibussowitsch if (view) PetscCall(VecView(z2, PETSC_VIEWER_DRAW_WORLD));
1189566063dSJacob Faibussowitsch PetscCall(VecAXPY(z1, -1.0, z2));
1199566063dSJacob Faibussowitsch PetscCall(VecNorm(z1, NORM_1, &enorm));
12048a46eb9SPierre Jolivet if (enorm > 1.e-11) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Error norm of |z1 - z2| %g\n", (double)enorm));
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