xref: /petsc/src/mat/tests/ex153.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown static char help[] = "This program illustrates the use of PETSc-fftw interface for sequential real DFT\n";
2c4762a1bSJed Brown #include <petscmat.h>
3c4762a1bSJed Brown #include <fftw3-mpi.h>
4*9371c9d4SSatish Balay int main(int argc, char **args) {
5c4762a1bSJed Brown   PetscMPIInt rank, size;
6c4762a1bSJed Brown   PetscInt    N0 = 10, N1 = 10, N2 = 10, N3 = 10, N4 = 10, N = N0 * N1 * N2 * N3 * N4;
7c4762a1bSJed Brown   PetscRandom rdm;
8c4762a1bSJed Brown   PetscReal   enorm;
9c4762a1bSJed Brown   Vec         x, y, z, input, output;
10c4762a1bSJed Brown   Mat         A;
11c4762a1bSJed Brown   PetscInt    DIM, dim[5], vsize;
12c4762a1bSJed Brown   PetscReal   fac;
13c4762a1bSJed Brown 
14327415f7SBarry Smith   PetscFunctionBeginUser;
159566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
16c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
17c4762a1bSJed Brown   SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "This example requires real numbers");
18c4762a1bSJed Brown #endif
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
21c4762a1bSJed Brown 
22be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uni-processor example only");
239566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
249566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
259566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &input));
269566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(input, N, N));
279566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(input));
289566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(input, rdm));
299566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(input, &output));
30c4762a1bSJed Brown 
31*9371c9d4SSatish Balay   DIM    = 5;
32*9371c9d4SSatish Balay   dim[0] = N0;
33*9371c9d4SSatish Balay   dim[1] = N1;
34*9371c9d4SSatish Balay   dim[2] = N2;
35*9371c9d4SSatish Balay   dim[3] = N3;
36*9371c9d4SSatish Balay   dim[4] = N4;
379566063dSJacob Faibussowitsch   PetscCall(MatCreateFFT(PETSC_COMM_SELF, DIM, dim, MATFFTW, &A));
389566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &x, &y));
399566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &z, NULL));
40c4762a1bSJed Brown 
419566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &vsize));
42c4762a1bSJed Brown   printf("The vector size  of input from the main routine is %d\n", vsize);
43c4762a1bSJed Brown 
449566063dSJacob Faibussowitsch   PetscCall(VecGetSize(z, &vsize));
45c4762a1bSJed Brown   printf("The vector size of output from the main routine is %d\n", vsize);
46c4762a1bSJed Brown 
479566063dSJacob Faibussowitsch   PetscCall(InputTransformFFT(A, input, x));
48c4762a1bSJed Brown 
499566063dSJacob Faibussowitsch   PetscCall(MatMult(A, x, y));
509566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A, y, z));
51c4762a1bSJed Brown 
529566063dSJacob Faibussowitsch   PetscCall(OutputTransformFFT(A, z, output));
53c4762a1bSJed Brown   fac = 1.0 / (PetscReal)N;
549566063dSJacob Faibussowitsch   PetscCall(VecScale(output, fac));
55c4762a1bSJed Brown   /*
569566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(input));
579566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(input));
589566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(output));
599566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(output));
60c4762a1bSJed Brown 
619566063dSJacob Faibussowitsch   PetscCall(VecView(input,PETSC_VIEWER_STDOUT_WORLD));
629566063dSJacob Faibussowitsch   PetscCall(VecView(output,PETSC_VIEWER_STDOUT_WORLD));
63c4762a1bSJed Brown */
649566063dSJacob Faibussowitsch   PetscCall(VecAXPY(output, -1.0, input));
659566063dSJacob Faibussowitsch   PetscCall(VecNorm(output, NORM_1, &enorm));
66c4762a1bSJed Brown   /*  if (enorm > 1.e-14) { */
679566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Error norm of |x - z| %e\n", enorm));
68c4762a1bSJed Brown   /*      } */
69c4762a1bSJed Brown 
709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&output));
719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&input));
729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&z));
759566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
769566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
779566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
78b122ec5aSJacob Faibussowitsch   return 0;
79c4762a1bSJed Brown }
80