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