Actual source code: ex153.c

  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>
  6: PetscInt main(PetscInt argc,char **args)
  7: {
  8:   PetscErrorCode  ierr;
  9:   PetscMPIInt     rank,size;
 10:   PetscInt        N0=10,N1=10,N2=10,N3=10,N4=10,N=N0*N1*N2*N3*N4;
 11:   PetscRandom     rdm;
 12:   PetscReal       enorm;
 13:   Vec             x,y,z,input,output;
 14:   Mat             A;
 15:   PetscInt        DIM, dim[5],vsize;
 16:   PetscReal       fac;

 18:   PetscInitialize(&argc,&args,(char *)0,help);
 19:   MPI_Comm_size(PETSC_COMM_WORLD, &size);
 20:   MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
 21: 
 22:   if (size!=1)
 23:     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This is a uni-processor example only");
 24:   PetscRandomCreate(PETSC_COMM_SELF, &rdm);
 25:   PetscRandomSetFromOptions(rdm);
 26:   VecCreate(PETSC_COMM_SELF,&input);
 27:   VecSetSizes(input,N,N);
 28:   VecSetFromOptions(input);
 29:   VecSetRandom(input,rdm);
 30:   VecDuplicate(input,&output);
 31: 
 32:   DIM = 5; dim[0] = N0; dim[1] = N1; dim[2] = N2; dim[3] = N3; dim[4] = N4;
 33:   MatCreateFFT(PETSC_COMM_SELF,DIM,dim,MATFFTW,&A);
 34:   MatGetVecs(A,&x,&y);
 35:   MatGetVecs(A,&z,PETSC_NULL);

 37:   VecGetSize(x,&vsize);
 38:   printf("The vector size  of input from the main routine is %d\n",vsize);

 40:   VecGetSize(z,&vsize);
 41:   printf("The vector size of output from the main routine is %d\n",vsize);

 43:   InputTransformFFT(A,input,x);

 45:   MatMult(A,x,y);
 46:   MatMultTranspose(A,y,z);

 48:   OutputTransformFFT(A,z,output);
 49:   fac = 1.0/(PetscReal)N;
 50:   VecScale(output,fac);
 51: /*  
 52:   VecAssemblyBegin(input);
 53:   VecAssemblyEnd(input);
 54:   VecAssemblyBegin(output);
 55:   VecAssemblyEnd(output);

 57:   VecView(input,PETSC_VIEWER_STDOUT_WORLD);
 58:   VecView(output,PETSC_VIEWER_STDOUT_WORLD);
 59: */
 60:   VecAXPY(output,-1.0,input);
 61:   VecNorm(output,NORM_1,&enorm);
 62: //  if (enorm > 1.e-14){
 63:       PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %e\n",enorm);
 64: //      }

 66:   VecDestroy(&output);
 67:   VecDestroy(&input);
 68:   VecDestroy(&x);
 69:   VecDestroy(&y);
 70:   VecDestroy(&z);
 71:   MatDestroy(&A);
 72:   PetscRandomDestroy(&rdm);
 73:   PetscFinalize();
 74:   return 0;

 76: }