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; 14 15 PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 16 #if defined(PETSC_USE_COMPLEX) 17 SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers"); 18 #endif 19 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 20 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 21 22 PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE, "This is a uni-processor example only"); 23 PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm)); 24 PetscCall(PetscRandomSetFromOptions(rdm)); 25 PetscCall(VecCreate(PETSC_COMM_SELF,&input)); 26 PetscCall(VecSetSizes(input,N,N)); 27 PetscCall(VecSetFromOptions(input)); 28 PetscCall(VecSetRandom(input,rdm)); 29 PetscCall(VecDuplicate(input,&output)); 30 31 DIM = 5; dim[0] = N0; dim[1] = N1; dim[2] = N2; dim[3] = N3; dim[4] = N4; 32 PetscCall(MatCreateFFT(PETSC_COMM_SELF,DIM,dim,MATFFTW,&A)); 33 PetscCall(MatCreateVecs(A,&x,&y)); 34 PetscCall(MatCreateVecs(A,&z,NULL)); 35 36 PetscCall(VecGetSize(x,&vsize)); 37 printf("The vector size of input from the main routine is %d\n",vsize); 38 39 PetscCall(VecGetSize(z,&vsize)); 40 printf("The vector size of output from the main routine is %d\n",vsize); 41 42 PetscCall(InputTransformFFT(A,input,x)); 43 44 PetscCall(MatMult(A,x,y)); 45 PetscCall(MatMultTranspose(A,y,z)); 46 47 PetscCall(OutputTransformFFT(A,z,output)); 48 fac = 1.0/(PetscReal)N; 49 PetscCall(VecScale(output,fac)); 50 /* 51 PetscCall(VecAssemblyBegin(input)); 52 PetscCall(VecAssemblyEnd(input)); 53 PetscCall(VecAssemblyBegin(output)); 54 PetscCall(VecAssemblyEnd(output)); 55 56 PetscCall(VecView(input,PETSC_VIEWER_STDOUT_WORLD)); 57 PetscCall(VecView(output,PETSC_VIEWER_STDOUT_WORLD)); 58 */ 59 PetscCall(VecAXPY(output,-1.0,input)); 60 PetscCall(VecNorm(output,NORM_1,&enorm)); 61 /* if (enorm > 1.e-14) { */ 62 PetscCall(PetscPrintf(PETSC_COMM_SELF," Error norm of |x - z| %e\n",enorm)); 63 /* } */ 64 65 PetscCall(VecDestroy(&output)); 66 PetscCall(VecDestroy(&input)); 67 PetscCall(VecDestroy(&x)); 68 PetscCall(VecDestroy(&y)); 69 PetscCall(VecDestroy(&z)); 70 PetscCall(MatDestroy(&A)); 71 PetscCall(PetscRandomDestroy(&rdm)); 72 PetscCall(PetscFinalize()); 73 return 0; 74 75 } 76