Actual source code: ex155.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: /*extern PetscErrorCode MatCreateVecsFFT(Mat,Vec *,Vec *,Vec *);*/
  5: int main(int argc, char **args)
  6: {
  7:   PetscMPIInt rank, size;
  8:   PetscInt    N0 = 4096, N1 = 4096, N2 = 256, N3 = 10, N4 = 10, N = N0 * N1;
  9:   PetscRandom rdm;
 10:   PetscReal   enorm;
 11:   Vec         x, y, z, input, output;
 12:   Mat         A;
 13:   PetscInt    DIM, dim[5], vsize, row, col;
 14:   PetscReal   fac;

 17:   PetscInitialize(&argc, &args, (char *)0, help);
 18:   MPI_Comm_size(PETSC_COMM_WORLD, &size);
 19:   MPI_Comm_rank(PETSC_COMM_WORLD, &rank);

 21: #if defined(PETSC_USE_COMPLEX)
 22:   SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Example for Real DFT. Your current data type is complex!");
 23: #endif
 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:   VecSetRandom(input, rdm);
 30:   VecDuplicate(input, &output);

 32:   DIM    = 2;
 33:   dim[0] = N0;
 34:   dim[1] = N1;
 35:   dim[2] = N2;
 36:   dim[3] = N3;
 37:   dim[4] = N4;
 38:   MatCreateFFT(PETSC_COMM_WORLD, DIM, dim, MATFFTW, &A);
 39:   MatGetLocalSize(A, &row, &col);
 40:   printf("The Matrix size  is %d and %d from process %d\n", row, col, rank);
 41:   MatCreateVecsFFTW(A, &x, &y, &z);

 43:   VecGetSize(x, &vsize);

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

 48:   VecScatterPetscToFFTW(A, input, x);
 49:   /*VecDestroy(&input);*/
 50:   MatMult(A, x, y);
 51:   MatMultTranspose(A, y, z);
 52:   VecScatterFFTWToPetsc(A, z, output);
 53:   /*VecDestroy(&z);*/
 54:   fac = 1.0 / (PetscReal)N;
 55:   VecScale(output, fac);

 57:   VecAssemblyBegin(input);
 58:   VecAssemblyEnd(input);
 59:   VecAssemblyBegin(output);
 60:   VecAssemblyEnd(output);

 62:   /*  VecView(input,PETSC_VIEWER_STDOUT_WORLD);*/
 63:   /*  VecView(output,PETSC_VIEWER_STDOUT_WORLD);*/

 65:   VecAXPY(output, -1.0, input);
 66:   VecNorm(output, NORM_1, &enorm);
 67:   if (enorm > 1.e-14) PetscPrintf(PETSC_COMM_SELF, "  Error norm of |x - z| %e\n", enorm);

 69:   VecDestroy(&x);
 70:   VecDestroy(&y);
 71:   VecDestroy(&z);
 72:   VecDestroy(&output);
 73:   VecDestroy(&input);
 74:   MatDestroy(&A);
 75:   PetscRandomDestroy(&rdm);
 76:   PetscFinalize();
 77:   return 0;
 78: }