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