xref: /petsc/src/mat/tests/ex149.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, N = N0 * N1 * N2;
11c4762a1bSJed Brown   PetscRandom rdm;
12c4762a1bSJed Brown   PetscScalar a;
13c4762a1bSJed Brown   PetscReal   enorm;
14c4762a1bSJed Brown   Vec         x, y, z, input, output;
15c4762a1bSJed Brown   PetscBool   view = PETSC_FALSE, use_interface = PETSC_TRUE;
16c4762a1bSJed Brown   Mat         A;
17c4762a1bSJed Brown   PetscInt    DIM, dim[3], vsize;
18c4762a1bSJed Brown   PetscReal   fac;
19c4762a1bSJed Brown 
20327415f7SBarry Smith   PetscFunctionBeginUser;
21*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
220fdf79fbSJacob Faibussowitsch   PetscCheck(!PetscDefined(USE_COMPLEX), PETSC_COMM_WORLD, PETSC_ERR_SUP, "This example requires real numbers");
23c4762a1bSJed Brown 
249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
26c4762a1bSJed Brown 
279566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
289566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
29c4762a1bSJed Brown 
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(VecDuplicate(input, &output));
359566063dSJacob Faibussowitsch   /*  PetscCall(VecGetSize(input,&vsize)); */
36c4762a1bSJed Brown   /*  printf("Size of the input Vector is %d\n",vsize); */
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   DIM    = 3;
399371c9d4SSatish Balay   dim[0] = N0;
409371c9d4SSatish Balay   dim[1] = N1;
419371c9d4SSatish Balay   dim[2] = N2;
42c4762a1bSJed Brown 
439566063dSJacob Faibussowitsch   PetscCall(MatCreateFFT(PETSC_COMM_WORLD, DIM, dim, MATFFTW, &A));
449566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &x, &y));
459566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &z, NULL));
469566063dSJacob Faibussowitsch   PetscCall(VecGetSize(y, &vsize));
47c4762a1bSJed Brown   printf("The vector size from the main routine is %d\n", vsize);
48c4762a1bSJed Brown 
499566063dSJacob Faibussowitsch   PetscCall(InputTransformFFT(A, input, x));
509566063dSJacob Faibussowitsch   PetscCall(MatMult(A, x, y));
519566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A, y, z));
529566063dSJacob Faibussowitsch   PetscCall(OutputTransformFFT(A, z, output));
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   fac = 1.0 / (PetscReal)N;
559566063dSJacob Faibussowitsch   PetscCall(VecScale(output, fac));
56c4762a1bSJed Brown 
579566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(input));
589566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(input));
599566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(output));
609566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(output));
61c4762a1bSJed Brown 
629566063dSJacob Faibussowitsch   PetscCall(VecView(input, PETSC_VIEWER_STDOUT_WORLD));
639566063dSJacob Faibussowitsch   PetscCall(VecView(output, PETSC_VIEWER_STDOUT_WORLD));
64c4762a1bSJed Brown 
659566063dSJacob Faibussowitsch   PetscCall(VecAXPY(output, -1.0, input));
669566063dSJacob Faibussowitsch   PetscCall(VecNorm(output, NORM_1, &enorm));
67c4762a1bSJed Brown   /*  if (enorm > 1.e-14) { */
6848a46eb9SPierre Jolivet   if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Error norm of |x - z| %e\n", enorm));
69c4762a1bSJed Brown   /*      } */
70c4762a1bSJed Brown 
719566063dSJacob Faibussowitsch   /* PetscCall(MatCreateVecs(A,&z,NULL)); */
72c4762a1bSJed Brown   /*  printf("Vector size from ex148 %d\n",vsize); */
739566063dSJacob Faibussowitsch   /*  PetscCall(PetscObjectSetName((PetscObject) x, "Real space vector")); */
749566063dSJacob Faibussowitsch   /*      PetscCall(PetscObjectSetName((PetscObject) y, "Frequency space vector")); */
759566063dSJacob Faibussowitsch   /*      PetscCall(PetscObjectSetName((PetscObject) z, "Reconstructed vector")); */
76c4762a1bSJed Brown 
779566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
78b122ec5aSJacob Faibussowitsch   return 0;
79c4762a1bSJed Brown }
80