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