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; 15 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; dim[0] = N0; dim[1] = N1; dim[2] = N2; dim[3] = N3; dim[4] = N4; 38 PetscCall(MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A)); 39 PetscCall(MatCreateVecsFFTW(A,&x,&y,&z)); 40 /* PetscCall(MatCreateVecs(A,&x,&y)); */ 41 /* PetscCall(MatCreateVecs(A,&z,NULL)); */ 42 43 PetscCall(VecGetSize(x,&vsize)); 44 printf("The vector size of input from the main routine is %d\n",vsize); 45 46 PetscCall(VecGetSize(z,&vsize)); 47 printf("The vector size of output from the main routine is %d\n",vsize); 48 49 PetscCall(InputTransformFFT(A,input,x)); 50 51 PetscCall(MatMult(A,x,y)); 52 PetscCall(VecAssemblyBegin(y)); 53 PetscCall(VecAssemblyEnd(y)); 54 PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); 55 56 PetscCall(MatMultTranspose(A,y,z)); 57 58 PetscCall(OutputTransformFFT(A,z,output)); 59 fac = 1.0/(PetscReal)N; 60 PetscCall(VecScale(output,fac)); 61 62 PetscCall(VecAssemblyBegin(input)); 63 PetscCall(VecAssemblyEnd(input)); 64 PetscCall(VecAssemblyBegin(output)); 65 PetscCall(VecAssemblyEnd(output)); 66 67 /* PetscCall(VecView(input,PETSC_VIEWER_STDOUT_WORLD)); */ 68 /* PetscCall(VecView(output,PETSC_VIEWER_STDOUT_WORLD)); */ 69 70 PetscCall(VecAXPY(output,-1.0,input)); 71 PetscCall(VecNorm(output,NORM_1,&enorm)); 72 /* if (enorm > 1.e-14) { */ 73 PetscCall(PetscPrintf(PETSC_COMM_SELF," Error norm of |x - z| %e\n",enorm)); 74 /* } */ 75 76 PetscCall(VecDestroy(&output)); 77 PetscCall(VecDestroy(&input)); 78 PetscCall(VecDestroy(&x)); 79 PetscCall(VecDestroy(&y)); 80 PetscCall(VecDestroy(&z)); 81 PetscCall(MatDestroy(&A)); 82 PetscCall(PetscRandomDestroy(&rdm)); 83 PetscCall(PetscFinalize()); 84 return 0; 85 86 } 87