1c4762a1bSJed Brown static char help[]="This program illustrates the use of PETSc-fftw interface for real 2D DFT.\n\ 2c4762a1bSJed Brown See ~petsc/src/mat/tests/ex158.c for general cases. \n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown PetscMPIInt rank,size; 9c4762a1bSJed Brown PetscInt N0=50,N1=50,N=N0*N1; 10c4762a1bSJed Brown PetscRandom rdm; 11c4762a1bSJed Brown PetscReal enorm; 12c4762a1bSJed Brown Vec x,y,z,input,output; 13c4762a1bSJed Brown Mat A; 14c4762a1bSJed Brown PetscInt DIM,dim[2]; 15c4762a1bSJed Brown PetscReal fac; 16c4762a1bSJed Brown 17*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 18c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 19c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers"); 20c4762a1bSJed Brown #endif 21c4762a1bSJed Brown 22*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 23*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 24c4762a1bSJed Brown 25c4762a1bSJed Brown /* Create and set PETSc vectors 'input' and 'output' */ 26*9566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm)); 27*9566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 28c4762a1bSJed Brown 29*9566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&input)); 30*9566063dSJacob Faibussowitsch PetscCall(VecSetSizes(input,PETSC_DECIDE,N0*N1)); 31*9566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(input)); 32*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(input,rdm)); 33*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(input,&output)); 34*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)input, "Real space vector")); 35*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)output, "Reconstructed vector")); 36c4762a1bSJed Brown 37c4762a1bSJed Brown /* Get FFTW vectors 'x', 'y' and 'z' */ 38c4762a1bSJed Brown DIM = 2; 39c4762a1bSJed Brown dim[0] = N0; dim[1] = N1; 40*9566063dSJacob Faibussowitsch PetscCall(MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A)); 41*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecsFFTW(A,&x,&y,&z)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Scatter PETSc vector 'input' to FFTW vector 'x' */ 44*9566063dSJacob Faibussowitsch PetscCall(VecScatterPetscToFFTW(A,input,x)); 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* Apply forward FFT */ 47*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,y)); 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* Apply backward FFT */ 50*9566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A,y,z)); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Scatter FFTW vector 'z' to PETSc vector 'output' */ 53*9566063dSJacob Faibussowitsch PetscCall(VecScatterFFTWToPetsc(A,z,output)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Check accuracy */ 56c4762a1bSJed Brown fac = 1.0/(PetscReal)N; 57*9566063dSJacob Faibussowitsch PetscCall(VecScale(output,fac)); 58*9566063dSJacob Faibussowitsch PetscCall(VecAXPY(output,-1.0,input)); 59*9566063dSJacob Faibussowitsch PetscCall(VecNorm(output,NORM_1,&enorm)); 60dd400576SPatrick Sanan if (enorm > 1.e-11 && rank == 0) { 61*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," Error norm of |x - z| %e\n",enorm)); 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* Free spaces */ 65*9566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 66*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&input)); 67*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&output)); 68*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 69*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 70*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&z)); 71*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 72c4762a1bSJed Brown 73*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 74b122ec5aSJacob Faibussowitsch return 0; 75c4762a1bSJed Brown } 76c4762a1bSJed Brown 77c4762a1bSJed Brown /*TEST 78c4762a1bSJed Brown 79c4762a1bSJed Brown build: 80c4762a1bSJed Brown requires: fftw !complex 81c4762a1bSJed Brown 82c4762a1bSJed Brown test: 83c4762a1bSJed Brown output_file: output/ex148.out 84c4762a1bSJed Brown 85c4762a1bSJed Brown test: 86c4762a1bSJed Brown suffix: 2 87c4762a1bSJed Brown nsize: 3 88c4762a1bSJed Brown output_file: output/ex148.out 89c4762a1bSJed Brown 90c4762a1bSJed Brown TEST*/ 91