1c4762a1bSJed Brown static char help[]="This program illustrates the use of PETSc-fftw interface for sequential real DFT\n"; 2c4762a1bSJed Brown #include <petscmat.h> 3c4762a1bSJed Brown #include <fftw3-mpi.h> 4ad0e0ea3SStefano Zampini int main(int argc,char **args) 5c4762a1bSJed Brown { 6c4762a1bSJed Brown PetscMPIInt rank,size; 7c4762a1bSJed Brown PetscInt N0=10,N1=10,N2=10,N3=10,N4=10,N=N0*N1*N2*N3*N4; 8c4762a1bSJed Brown PetscRandom rdm; 9c4762a1bSJed Brown PetscReal enorm; 10c4762a1bSJed Brown Vec x,y,z,input,output; 11c4762a1bSJed Brown Mat A; 12c4762a1bSJed Brown PetscInt DIM, dim[5],vsize; 13c4762a1bSJed Brown PetscReal fac; 14c4762a1bSJed Brown 15*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 16c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 17c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers"); 18c4762a1bSJed Brown #endif 19*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 20*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 21c4762a1bSJed Brown 222c71b3e2SJacob Faibussowitsch PetscCheckFalse(size!=1,PETSC_COMM_WORLD,PETSC_ERR_SUP, "This is a uni-processor example only"); 23*9566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm)); 24*9566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 25*9566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF,&input)); 26*9566063dSJacob Faibussowitsch PetscCall(VecSetSizes(input,N,N)); 27*9566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(input)); 28*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(input,rdm)); 29*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(input,&output)); 30c4762a1bSJed Brown 31c4762a1bSJed Brown DIM = 5; dim[0] = N0; dim[1] = N1; dim[2] = N2; dim[3] = N3; dim[4] = N4; 32*9566063dSJacob Faibussowitsch PetscCall(MatCreateFFT(PETSC_COMM_SELF,DIM,dim,MATFFTW,&A)); 33*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&x,&y)); 34*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&z,NULL)); 35c4762a1bSJed Brown 36*9566063dSJacob Faibussowitsch PetscCall(VecGetSize(x,&vsize)); 37c4762a1bSJed Brown printf("The vector size of input from the main routine is %d\n",vsize); 38c4762a1bSJed Brown 39*9566063dSJacob Faibussowitsch PetscCall(VecGetSize(z,&vsize)); 40c4762a1bSJed Brown printf("The vector size of output from the main routine is %d\n",vsize); 41c4762a1bSJed Brown 42*9566063dSJacob Faibussowitsch PetscCall(InputTransformFFT(A,input,x)); 43c4762a1bSJed Brown 44*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,y)); 45*9566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A,y,z)); 46c4762a1bSJed Brown 47*9566063dSJacob Faibussowitsch PetscCall(OutputTransformFFT(A,z,output)); 48c4762a1bSJed Brown fac = 1.0/(PetscReal)N; 49*9566063dSJacob Faibussowitsch PetscCall(VecScale(output,fac)); 50c4762a1bSJed Brown /* 51*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(input)); 52*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(input)); 53*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(output)); 54*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(output)); 55c4762a1bSJed Brown 56*9566063dSJacob Faibussowitsch PetscCall(VecView(input,PETSC_VIEWER_STDOUT_WORLD)); 57*9566063dSJacob Faibussowitsch PetscCall(VecView(output,PETSC_VIEWER_STDOUT_WORLD)); 58c4762a1bSJed Brown */ 59*9566063dSJacob Faibussowitsch PetscCall(VecAXPY(output,-1.0,input)); 60*9566063dSJacob Faibussowitsch PetscCall(VecNorm(output,NORM_1,&enorm)); 61c4762a1bSJed Brown /* if (enorm > 1.e-14) { */ 62*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," Error norm of |x - z| %e\n",enorm)); 63c4762a1bSJed Brown /* } */ 64c4762a1bSJed Brown 65*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&output)); 66*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&input)); 67*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 68*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 69*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&z)); 70*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 71*9566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 72*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 73b122ec5aSJacob Faibussowitsch return 0; 74c4762a1bSJed Brown 75c4762a1bSJed Brown } 76