xref: /petsc/src/mat/tests/ex155.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
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>
4c4762a1bSJed Brown /*extern PetscErrorCode MatCreateVecsFFT(Mat,Vec *,Vec *,Vec *);*/
main(int argc,char ** args)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown   PetscMPIInt rank, size;
8c4762a1bSJed Brown   PetscInt    N0 = 4096, N1 = 4096, N2 = 256, N3 = 10, N4 = 10, N = N0 * N1;
9c4762a1bSJed Brown   PetscRandom rdm;
10c4762a1bSJed Brown   PetscReal   enorm;
11c4762a1bSJed Brown   Vec         x, y, z, input, output;
12c4762a1bSJed Brown   Mat         A;
13c4762a1bSJed Brown   PetscInt    DIM, dim[5], vsize, row, col;
14c4762a1bSJed Brown   PetscReal   fac;
15c4762a1bSJed Brown 
16327415f7SBarry Smith   PetscFunctionBeginUser;
17*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
20c4762a1bSJed Brown 
21c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
22c4762a1bSJed Brown   SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Example for Real DFT. Your current data type is complex!");
23c4762a1bSJed Brown #endif
249566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
259566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
269566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &input));
279566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(input, PETSC_DECIDE, N));
289566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(input));
299566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(input, rdm));
309566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(input, &output));
31c4762a1bSJed Brown 
329371c9d4SSatish Balay   DIM    = 2;
339371c9d4SSatish Balay   dim[0] = N0;
349371c9d4SSatish Balay   dim[1] = N1;
359371c9d4SSatish Balay   dim[2] = N2;
369371c9d4SSatish Balay   dim[3] = N3;
379371c9d4SSatish Balay   dim[4] = N4;
389566063dSJacob Faibussowitsch   PetscCall(MatCreateFFT(PETSC_COMM_WORLD, DIM, dim, MATFFTW, &A));
399566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &row, &col));
40c4762a1bSJed Brown   printf("The Matrix size  is %d and %d from process %d\n", row, col, rank);
419566063dSJacob Faibussowitsch   PetscCall(MatCreateVecsFFTW(A, &x, &y, &z));
42c4762a1bSJed Brown 
439566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &vsize));
44c4762a1bSJed Brown 
459566063dSJacob Faibussowitsch   PetscCall(VecGetSize(z, &vsize));
46c4762a1bSJed Brown   printf("The vector size of output from the main routine is %d\n", vsize);
47c4762a1bSJed Brown 
489566063dSJacob Faibussowitsch   PetscCall(VecScatterPetscToFFTW(A, input, x));
499566063dSJacob Faibussowitsch   /*PetscCall(VecDestroy(&input));*/
509566063dSJacob Faibussowitsch   PetscCall(MatMult(A, x, y));
519566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A, y, z));
529566063dSJacob Faibussowitsch   PetscCall(VecScatterFFTWToPetsc(A, z, output));
539566063dSJacob Faibussowitsch   /*PetscCall(VecDestroy(&z));*/
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));
6748a46eb9SPierre Jolivet   if (enorm > 1.e-14) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Error norm of |x - z| %e\n", enorm));
68c4762a1bSJed Brown 
699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&z));
729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&output));
739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&input));
749566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
759566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
769566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
77b122ec5aSJacob Faibussowitsch   return 0;
78c4762a1bSJed Brown }
79