xref: /petsc/src/mat/tests/ex155.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 static char help[]="This program illustrates the use of PETSc-fftw interface for parallel real DFT\n";
2 #include <petscmat.h>
3 #include <fftw3-mpi.h>
4 /*extern PetscErrorCode MatCreateVecsFFT(Mat,Vec *,Vec *,Vec *);*/
5 int main(int argc,char **args)
6 {
7   PetscMPIInt    rank,size;
8   PetscInt       N0=4096,N1=4096,N2=256,N3=10,N4=10,N=N0*N1;
9   PetscRandom    rdm;
10   PetscReal      enorm;
11   Vec            x,y,z,input,output;
12   Mat            A;
13   PetscInt       DIM, dim[5],vsize,row,col;
14   PetscReal      fac;
15 
16   PetscFunctionBeginUser;
17   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
18   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
19   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
20 
21 #if defined(PETSC_USE_COMPLEX)
22   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "Example for Real DFT. Your current data type is complex!");
23 #endif
24   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
25   PetscCall(PetscRandomSetFromOptions(rdm));
26   PetscCall(VecCreate(PETSC_COMM_WORLD,&input));
27   PetscCall(VecSetSizes(input,PETSC_DECIDE,N));
28   PetscCall(VecSetFromOptions(input));
29   PetscCall(VecSetRandom(input,rdm));
30   PetscCall(VecDuplicate(input,&output));
31 
32   DIM  = 2; dim[0] = N0; dim[1] = N1; dim[2] = N2; dim[3] = N3; dim[4] = N4;
33   PetscCall(MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A));
34   PetscCall(MatGetLocalSize(A,&row,&col));
35   printf("The Matrix size  is %d and %d from process %d\n",row,col,rank);
36   PetscCall(MatCreateVecsFFTW(A,&x,&y,&z));
37 
38   PetscCall(VecGetSize(x,&vsize));
39 
40   PetscCall(VecGetSize(z,&vsize));
41   printf("The vector size of output from the main routine is %d\n",vsize);
42 
43   PetscCall(VecScatterPetscToFFTW(A,input,x));
44   /*PetscCall(VecDestroy(&input));*/
45   PetscCall(MatMult(A,x,y));
46   PetscCall(MatMultTranspose(A,y,z));
47   PetscCall(VecScatterFFTWToPetsc(A,z,output));
48   /*PetscCall(VecDestroy(&z));*/
49   fac  = 1.0/(PetscReal)N;
50   PetscCall(VecScale(output,fac));
51 
52   PetscCall(VecAssemblyBegin(input));
53   PetscCall(VecAssemblyEnd(input));
54   PetscCall(VecAssemblyBegin(output));
55   PetscCall(VecAssemblyEnd(output));
56 
57 /*  PetscCall(VecView(input,PETSC_VIEWER_STDOUT_WORLD));*/
58 /*  PetscCall(VecView(output,PETSC_VIEWER_STDOUT_WORLD));*/
59 
60   PetscCall(VecAXPY(output,-1.0,input));
61   PetscCall(VecNorm(output,NORM_1,&enorm));
62   if (enorm > 1.e-14) {
63     PetscCall(PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %e\n",enorm));
64   }
65 
66   PetscCall(VecDestroy(&x));
67   PetscCall(VecDestroy(&y));
68   PetscCall(VecDestroy(&z));
69   PetscCall(VecDestroy(&output));
70   PetscCall(VecDestroy(&input));
71   PetscCall(MatDestroy(&A));
72   PetscCall(PetscRandomDestroy(&rdm));
73   PetscCall(PetscFinalize());
74   return 0;
75 
76 }
77