xref: /petsc/src/mat/tests/ex157.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown static char help[] = "This program illustrates the use of PETSc-fftw interface for parallel 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 = 2048, N1 = 2048, N2 = 3, N3 = 5, N4 = 5, N = N0 * N1;
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   PetscScalar one = 1, two = 2, three = 3;
14c4762a1bSJed Brown 
15327415f7SBarry Smith   PetscFunctionBeginUser;
169566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
17c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
18c4762a1bSJed Brown   SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "This example requires real numbers");
19c4762a1bSJed Brown #endif
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
22c4762a1bSJed Brown 
239566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
249566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
259566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &input));
269566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(input, PETSC_DECIDE, N));
279566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(input));
289566063dSJacob Faibussowitsch   /*  PetscCall(VecSet(input,one)); */
299566063dSJacob Faibussowitsch   /*  PetscCall(VecSetValue(input,1,two,INSERT_VALUES)); */
309566063dSJacob Faibussowitsch   /*  PetscCall(VecSetValue(input,2,three,INSERT_VALUES)); */
319566063dSJacob Faibussowitsch   /*  PetscCall(VecSetValue(input,3,three,INSERT_VALUES)); */
329566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(input, rdm));
339566063dSJacob Faibussowitsch   /*  PetscCall(VecSetRandom(input,rdm)); */
349566063dSJacob Faibussowitsch   /*  PetscCall(VecSetRandom(input,rdm)); */
359566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(input, &output));
36c4762a1bSJed Brown 
37*9371c9d4SSatish Balay   DIM    = 2;
38*9371c9d4SSatish Balay   dim[0] = N0;
39*9371c9d4SSatish Balay   dim[1] = N1;
40*9371c9d4SSatish Balay   dim[2] = N2;
41*9371c9d4SSatish Balay   dim[3] = N3;
42*9371c9d4SSatish Balay   dim[4] = N4;
439566063dSJacob Faibussowitsch   PetscCall(MatCreateFFT(PETSC_COMM_WORLD, DIM, dim, MATFFTW, &A));
449566063dSJacob Faibussowitsch   PetscCall(MatCreateVecsFFTW(A, &x, &y, &z));
459566063dSJacob Faibussowitsch   /*  PetscCall(MatCreateVecs(A,&x,&y)); */
469566063dSJacob Faibussowitsch   /*  PetscCall(MatCreateVecs(A,&z,NULL)); */
47c4762a1bSJed Brown 
489566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &vsize));
49c4762a1bSJed Brown   printf("The vector size  of input from the main routine is %d\n", vsize);
50c4762a1bSJed Brown 
519566063dSJacob Faibussowitsch   PetscCall(VecGetSize(z, &vsize));
52c4762a1bSJed Brown   printf("The vector size of output from the main routine is %d\n", vsize);
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(InputTransformFFT(A, input, x));
55c4762a1bSJed Brown 
569566063dSJacob Faibussowitsch   PetscCall(MatMult(A, x, y));
579566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(y));
589566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(y));
599566063dSJacob Faibussowitsch   PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
60c4762a1bSJed Brown 
619566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A, y, z));
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch   PetscCall(OutputTransformFFT(A, z, output));
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) { */
789566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Error norm of |x - z| %e\n", enorm));
79c4762a1bSJed Brown   /*      } */
80c4762a1bSJed Brown 
819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&output));
829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&input));
839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
859566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&z));
869566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
879566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
889566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
89b122ec5aSJacob Faibussowitsch   return 0;
90c4762a1bSJed Brown }
91