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>
  4: int main(int argc, char **args)
  5: {
  6:   PetscMPIInt rank, size;
  7:   PetscInt    N0 = 10, N1 = 10, N2 = 10, N3 = 10, N4 = 10, N = N0 * N1 * N2 * N3 * N4;
  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;

 15:   PetscFunctionBeginUser;
 16:   PetscCall(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:   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
 21:   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));

 23:   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uni-processor example only");
 24:   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
 25:   PetscCall(PetscRandomSetFromOptions(rdm));
 26:   PetscCall(VecCreate(PETSC_COMM_SELF, &input));
 27:   PetscCall(VecSetSizes(input, N, N));
 28:   PetscCall(VecSetFromOptions(input));
 29:   PetscCall(VecSetRandom(input, rdm));
 30:   PetscCall(VecDuplicate(input, &output));

 32:   DIM    = 5;
 33:   dim[0] = N0;
 34:   dim[1] = N1;
 35:   dim[2] = N2;
 36:   dim[3] = N3;
 37:   dim[4] = N4;
 38:   PetscCall(MatCreateFFT(PETSC_COMM_SELF, DIM, dim, MATFFTW, &A));
 39:   PetscCall(MatCreateVecs(A, &x, &y));
 40:   PetscCall(MatCreateVecs(A, &z, NULL));

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

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

 48:   PetscCall(InputTransformFFT(A, input, x));

 50:   PetscCall(MatMult(A, x, y));
 51:   PetscCall(MatMultTranspose(A, y, z));

 53:   PetscCall(OutputTransformFFT(A, z, output));
 54:   fac = 1.0 / (PetscReal)N;
 55:   PetscCall(VecScale(output, fac));
 56:   /*
 57:   PetscCall(VecAssemblyBegin(input));
 58:   PetscCall(VecAssemblyEnd(input));
 59:   PetscCall(VecAssemblyBegin(output));
 60:   PetscCall(VecAssemblyEnd(output));

 62:   PetscCall(VecView(input,PETSC_VIEWER_STDOUT_WORLD));
 63:   PetscCall(VecView(output,PETSC_VIEWER_STDOUT_WORLD));
 64: */
 65:   PetscCall(VecAXPY(output, -1.0, input));
 66:   PetscCall(VecNorm(output, NORM_1, &enorm));
 67:   /*  if (enorm > 1.e-14) { */
 68:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Error norm of |x - z| %e\n", enorm));
 69:   /*      } */

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