1 static char help[] = "This program illustrates the use of PETSc-fftw interface for parallel real DFT\n";
2 #include <petscmat.h>
3 #include <fftw3-mpi.h>
main(int argc,char ** args)4 int main(int argc, char **args)
5 {
6 PetscMPIInt rank, size;
7 PetscInt N0 = 2048, N1 = 2048, N2 = 3, N3 = 5, N4 = 5, N = N0 * N1;
8 PetscRandom rdm;
9 PetscReal enorm;
10 Vec x, y, z, input, output;
11 Mat A;
12 PetscInt DIM, dim[5], vsize;
13 PetscReal fac;
14 PetscScalar one = 1, two = 2, three = 3;
15
16 PetscFunctionBeginUser;
17 PetscCall(PetscInitialize(&argc, &args, NULL, help));
18 #if defined(PETSC_USE_COMPLEX)
19 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "This example requires real numbers");
20 #endif
21 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
22 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
23
24 PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
25 PetscCall(PetscRandomSetFromOptions(rdm));
26 PetscCall(VecCreate(PETSC_COMM_WORLD, &input));
27 PetscCall(VecSetSizes(input, PETSC_DECIDE, N));
28 PetscCall(VecSetFromOptions(input));
29 /* PetscCall(VecSet(input,one)); */
30 /* PetscCall(VecSetValue(input,1,two,INSERT_VALUES)); */
31 /* PetscCall(VecSetValue(input,2,three,INSERT_VALUES)); */
32 /* PetscCall(VecSetValue(input,3,three,INSERT_VALUES)); */
33 PetscCall(VecSetRandom(input, rdm));
34 /* PetscCall(VecSetRandom(input,rdm)); */
35 /* PetscCall(VecSetRandom(input,rdm)); */
36 PetscCall(VecDuplicate(input, &output));
37
38 DIM = 2;
39 dim[0] = N0;
40 dim[1] = N1;
41 dim[2] = N2;
42 dim[3] = N3;
43 dim[4] = N4;
44 PetscCall(MatCreateFFT(PETSC_COMM_WORLD, DIM, dim, MATFFTW, &A));
45 PetscCall(MatCreateVecsFFTW(A, &x, &y, &z));
46 /* PetscCall(MatCreateVecs(A,&x,&y)); */
47 /* PetscCall(MatCreateVecs(A,&z,NULL)); */
48
49 PetscCall(VecGetSize(x, &vsize));
50 printf("The vector size of input from the main routine is %d\n", vsize);
51
52 PetscCall(VecGetSize(z, &vsize));
53 printf("The vector size of output from the main routine is %d\n", vsize);
54
55 PetscCall(InputTransformFFT(A, input, x));
56
57 PetscCall(MatMult(A, x, y));
58 PetscCall(VecAssemblyBegin(y));
59 PetscCall(VecAssemblyEnd(y));
60 PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
61
62 PetscCall(MatMultTranspose(A, y, z));
63
64 PetscCall(OutputTransformFFT(A, z, output));
65 fac = 1.0 / (PetscReal)N;
66 PetscCall(VecScale(output, fac));
67
68 PetscCall(VecAssemblyBegin(input));
69 PetscCall(VecAssemblyEnd(input));
70 PetscCall(VecAssemblyBegin(output));
71 PetscCall(VecAssemblyEnd(output));
72
73 /* PetscCall(VecView(input,PETSC_VIEWER_STDOUT_WORLD)); */
74 /* PetscCall(VecView(output,PETSC_VIEWER_STDOUT_WORLD)); */
75
76 PetscCall(VecAXPY(output, -1.0, input));
77 PetscCall(VecNorm(output, NORM_1, &enorm));
78 /* if (enorm > 1.e-14) { */
79 PetscCall(PetscPrintf(PETSC_COMM_SELF, " Error norm of |x - z| %e\n", enorm));
80 /* } */
81
82 PetscCall(VecDestroy(&output));
83 PetscCall(VecDestroy(&input));
84 PetscCall(VecDestroy(&x));
85 PetscCall(VecDestroy(&y));
86 PetscCall(VecDestroy(&z));
87 PetscCall(MatDestroy(&A));
88 PetscCall(PetscRandomDestroy(&rdm));
89 PetscCall(PetscFinalize());
90 return 0;
91 }
92