xref: /petsc/src/mat/tests/ex150.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1c4762a1bSJed Brown static char help[] = "This program illustrates the use of PETSc-fftw interface for real DFT\n";
2c4762a1bSJed Brown #include <petscmat.h>
3c4762a1bSJed Brown #include <fftw3-mpi.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown extern PetscErrorCode InputTransformFFT(Mat, Vec, Vec);
6c4762a1bSJed Brown extern PetscErrorCode OutputTransformFFT(Mat, Vec, Vec);
main(int argc,char ** args)7d71ae5a4SJacob Faibussowitsch int                   main(int argc, char **args)
8d71ae5a4SJacob Faibussowitsch {
9c4762a1bSJed Brown   PetscMPIInt rank, size;
10c4762a1bSJed Brown   PetscInt    N0 = 3, N1 = 3, N2 = 3, N3 = 3, N4 = 3, N = N0 * N1 * N2 * N3;
11c4762a1bSJed Brown   PetscRandom rdm;
12c4762a1bSJed Brown   PetscReal   enorm;
13c4762a1bSJed Brown   Vec         x, y, z, input, output;
14c4762a1bSJed Brown   Mat         A;
15c4762a1bSJed Brown   PetscInt    DIM, dim[5], vsize;
16c4762a1bSJed Brown   PetscReal   fac;
17c4762a1bSJed Brown   PetscScalar one = 1;
18c4762a1bSJed Brown 
19327415f7SBarry Smith   PetscFunctionBeginUser;
20*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
21c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
22c4762a1bSJed Brown   SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "This example requires real numbers");
23c4762a1bSJed Brown #endif
24c4762a1bSJed Brown 
259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
27c4762a1bSJed Brown 
289566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
299566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
309566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &input));
319566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(input, PETSC_DECIDE, N));
329566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(input));
339566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(input, rdm));
349566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(input));
359566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(input));
369566063dSJacob Faibussowitsch   /*  PetscCall(VecView(input,PETSC_VIEWER_STDOUT_WORLD)); */
379566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(input, &output));
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   DIM    = 4;
409371c9d4SSatish Balay   dim[0] = N0;
419371c9d4SSatish Balay   dim[1] = N1;
429371c9d4SSatish Balay   dim[2] = N2;
439371c9d4SSatish Balay   dim[3] = N3;
449371c9d4SSatish Balay   dim[4] = N4;
45c4762a1bSJed Brown 
469566063dSJacob Faibussowitsch   PetscCall(MatCreateFFT(PETSC_COMM_WORLD, DIM, dim, MATFFTW, &A));
479566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &x, &y));
489566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &z, NULL));
499566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &vsize));
50c4762a1bSJed Brown   printf("The vector size from the main routine is %d\n", vsize);
51c4762a1bSJed Brown 
529566063dSJacob Faibussowitsch   PetscCall(InputTransformFFT(A, input, x));
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(MatMult(A, x, y));
559566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(y));
569566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(y));
579566063dSJacob Faibussowitsch   PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
589566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A, y, z));
599566063dSJacob Faibussowitsch   PetscCall(VecGetSize(z, &vsize));
60c4762a1bSJed Brown   printf("The vector size of zfrom the main routine is %d\n", vsize);
61c4762a1bSJed Brown 
629566063dSJacob Faibussowitsch   PetscCall(OutputTransformFFT(A, z, output));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   fac = 1.0 / (PetscReal)N;
659566063dSJacob Faibussowitsch   PetscCall(VecScale(output, fac));
66c4762a1bSJed Brown 
679566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(input));
689566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(input));
699566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(output));
709566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(output));
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch   PetscCall(VecView(input, PETSC_VIEWER_STDOUT_WORLD));
739566063dSJacob Faibussowitsch   PetscCall(VecView(output, PETSC_VIEWER_STDOUT_WORLD));
74c4762a1bSJed Brown 
759566063dSJacob Faibussowitsch   PetscCall(VecAXPY(output, -1.0, input));
769566063dSJacob Faibussowitsch   PetscCall(VecNorm(output, NORM_1, &enorm));
77c4762a1bSJed Brown   /*  if (enorm > 1.e-14) { */
7848a46eb9SPierre Jolivet   if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Error norm of |x - z| %e\n", enorm));
79c4762a1bSJed Brown   /*      } */
80c4762a1bSJed Brown 
819566063dSJacob Faibussowitsch   /* PetscCall(MatCreateVecs(A,&z,NULL)); */
82c4762a1bSJed Brown   /*  printf("Vector size from ex148 %d\n",vsize); */
839566063dSJacob Faibussowitsch   /*  PetscCall(PetscObjectSetName((PetscObject) x, "Real space vector")); */
849566063dSJacob Faibussowitsch   /*      PetscCall(PetscObjectSetName((PetscObject) y, "Frequency space vector")); */
859566063dSJacob Faibussowitsch   /*      PetscCall(PetscObjectSetName((PetscObject) z, "Reconstructed vector")); */
86c4762a1bSJed Brown 
879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&output));
889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&input));
899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&z));
929566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
939566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
949566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
95b122ec5aSJacob Faibussowitsch   return 0;
96c4762a1bSJed Brown }
97