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;

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

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

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

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

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

 55:   InputTransformFFT(A, input, x);

 57:   MatMult(A, x, y);
 58:   VecAssemblyBegin(y);
 59:   VecAssemblyEnd(y);
 60:   VecView(y, PETSC_VIEWER_STDOUT_WORLD);

 62:   MatMultTranspose(A, y, z);

 64:   OutputTransformFFT(A, z, output);
 65:   fac = 1.0 / (PetscReal)N;
 66:   VecScale(output, fac);

 68:   VecAssemblyBegin(input);
 69:   VecAssemblyEnd(input);
 70:   VecAssemblyBegin(output);
 71:   VecAssemblyEnd(output);

 73:   /*  VecView(input,PETSC_VIEWER_STDOUT_WORLD); */
 74:   /*  VecView(output,PETSC_VIEWER_STDOUT_WORLD); */

 76:   VecAXPY(output, -1.0, input);
 77:   VecNorm(output, NORM_1, &enorm);
 78:   /*  if (enorm > 1.e-14) { */
 79:   PetscPrintf(PETSC_COMM_SELF, "  Error norm of |x - z| %e\n", enorm);
 80:   /*      } */

 82:   VecDestroy(&output);
 83:   VecDestroy(&input);
 84:   VecDestroy(&x);
 85:   VecDestroy(&y);
 86:   VecDestroy(&z);
 87:   MatDestroy(&A);
 88:   PetscRandomDestroy(&rdm);
 89:   PetscFinalize();
 90:   return 0;
 91: }