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 PetscMPIInt rank, size; 6 PetscInt N0 = 2048, N1 = 2048, N2 = 3, N3 = 5, N4 = 5, N = N0 * N1; 7 PetscRandom rdm; 8 PetscReal enorm; 9 Vec x, y, z, input, output; 10 Mat A; 11 PetscInt DIM, dim[5], vsize; 12 PetscReal fac; 13 PetscScalar one = 1, two = 2, three = 3; 14 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)); 22 23 PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm)); 24 PetscCall(PetscRandomSetFromOptions(rdm)); 25 PetscCall(VecCreate(PETSC_COMM_WORLD, &input)); 26 PetscCall(VecSetSizes(input, PETSC_DECIDE, N)); 27 PetscCall(VecSetFromOptions(input)); 28 /* PetscCall(VecSet(input,one)); */ 29 /* PetscCall(VecSetValue(input,1,two,INSERT_VALUES)); */ 30 /* PetscCall(VecSetValue(input,2,three,INSERT_VALUES)); */ 31 /* PetscCall(VecSetValue(input,3,three,INSERT_VALUES)); */ 32 PetscCall(VecSetRandom(input, rdm)); 33 /* PetscCall(VecSetRandom(input,rdm)); */ 34 /* PetscCall(VecSetRandom(input,rdm)); */ 35 PetscCall(VecDuplicate(input, &output)); 36 37 DIM = 2; 38 dim[0] = N0; 39 dim[1] = N1; 40 dim[2] = N2; 41 dim[3] = N3; 42 dim[4] = N4; 43 PetscCall(MatCreateFFT(PETSC_COMM_WORLD, DIM, dim, MATFFTW, &A)); 44 PetscCall(MatCreateVecsFFTW(A, &x, &y, &z)); 45 /* PetscCall(MatCreateVecs(A,&x,&y)); */ 46 /* PetscCall(MatCreateVecs(A,&z,NULL)); */ 47 48 PetscCall(VecGetSize(x, &vsize)); 49 printf("The vector size of input from the main routine is %d\n", vsize); 50 51 PetscCall(VecGetSize(z, &vsize)); 52 printf("The vector size of output from the main routine is %d\n", vsize); 53 54 PetscCall(InputTransformFFT(A, input, x)); 55 56 PetscCall(MatMult(A, x, y)); 57 PetscCall(VecAssemblyBegin(y)); 58 PetscCall(VecAssemblyEnd(y)); 59 PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD)); 60 61 PetscCall(MatMultTranspose(A, y, z)); 62 63 PetscCall(OutputTransformFFT(A, z, output)); 64 fac = 1.0 / (PetscReal)N; 65 PetscCall(VecScale(output, fac)); 66 67 PetscCall(VecAssemblyBegin(input)); 68 PetscCall(VecAssemblyEnd(input)); 69 PetscCall(VecAssemblyBegin(output)); 70 PetscCall(VecAssemblyEnd(output)); 71 72 /* PetscCall(VecView(input,PETSC_VIEWER_STDOUT_WORLD)); */ 73 /* PetscCall(VecView(output,PETSC_VIEWER_STDOUT_WORLD)); */ 74 75 PetscCall(VecAXPY(output, -1.0, input)); 76 PetscCall(VecNorm(output, NORM_1, &enorm)); 77 /* if (enorm > 1.e-14) { */ 78 PetscCall(PetscPrintf(PETSC_COMM_SELF, " Error norm of |x - z| %e\n", enorm)); 79 /* } */ 80 81 PetscCall(VecDestroy(&output)); 82 PetscCall(VecDestroy(&input)); 83 PetscCall(VecDestroy(&x)); 84 PetscCall(VecDestroy(&y)); 85 PetscCall(VecDestroy(&z)); 86 PetscCall(MatDestroy(&A)); 87 PetscCall(PetscRandomDestroy(&rdm)); 88 PetscCall(PetscFinalize()); 89 return 0; 90 } 91