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