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