Actual source code: ex157.c

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

 16:   PetscInitialize(&argc,&args,(char*)0,help);
 17: #if defined(PETSC_USE_COMPLEX)
 18:   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers");
 19: #endif
 20:   MPI_Comm_size(PETSC_COMM_WORLD, &size);
 21:   MPI_Comm_rank(PETSC_COMM_WORLD, &rank);

 23:   PetscRandomCreate(PETSC_COMM_WORLD, &rdm);
 24:   PetscRandomSetFromOptions(rdm);
 25:   VecCreate(PETSC_COMM_WORLD,&input);
 26:   VecSetSizes(input,PETSC_DECIDE,N);
 27:   VecSetFromOptions(input);
 28: /*  VecSet(input,one); */
 29: /*  VecSetValue(input,1,two,INSERT_VALUES); */
 30: /*  VecSetValue(input,2,three,INSERT_VALUES); */
 31: /*  VecSetValue(input,3,three,INSERT_VALUES); */
 32:   VecSetRandom(input,rdm);
 33: /*  VecSetRandom(input,rdm); */
 34: /*  VecSetRandom(input,rdm); */
 35:   VecDuplicate(input,&output);

 37:   DIM  = 2; dim[0] = N0; dim[1] = N1; dim[2] = N2; dim[3] = N3; dim[4] = N4;
 38:   MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A);
 39:   MatCreateVecsFFTW(A,&x,&y,&z);
 40: /*  MatCreateVecs(A,&x,&y); */
 41: /*  MatCreateVecs(A,&z,NULL); */

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

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

 49:   InputTransformFFT(A,input,x);

 51:   MatMult(A,x,y);
 52:   VecAssemblyBegin(y);
 53:   VecAssemblyEnd(y);
 54:   VecView(y,PETSC_VIEWER_STDOUT_WORLD);

 56:   MatMultTranspose(A,y,z);

 58:   OutputTransformFFT(A,z,output);
 59:   fac  = 1.0/(PetscReal)N;
 60:   VecScale(output,fac);

 62:   VecAssemblyBegin(input);
 63:   VecAssemblyEnd(input);
 64:   VecAssemblyBegin(output);
 65:   VecAssemblyEnd(output);

 67: /*  VecView(input,PETSC_VIEWER_STDOUT_WORLD); */
 68: /*  VecView(output,PETSC_VIEWER_STDOUT_WORLD); */

 70:   VecAXPY(output,-1.0,input);
 71:   VecNorm(output,NORM_1,&enorm);
 72: /*  if (enorm > 1.e-14) { */
 73:   PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %e\n",enorm);
 74: /*      } */

 76:   VecDestroy(&output);
 77:   VecDestroy(&input);
 78:   VecDestroy(&x);
 79:   VecDestroy(&y);
 80:   VecDestroy(&z);
 81:   MatDestroy(&A);
 82:   PetscRandomDestroy(&rdm);
 83:   PetscFinalize();
 84:   return 0;

 86: }