1c4762a1bSJed Brown static char help[]="This program illustrates the use of PETSc-fftw interface for real DFT\n"; 2c4762a1bSJed Brown #include <petscmat.h> 3c4762a1bSJed Brown #include <fftw3-mpi.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown extern PetscErrorCode InputTransformFFT(Mat,Vec,Vec); 6c4762a1bSJed Brown extern PetscErrorCode OutputTransformFFT(Mat,Vec,Vec); 7ad0e0ea3SStefano Zampini int main(int argc,char **args) 8c4762a1bSJed Brown { 9c4762a1bSJed Brown PetscErrorCode ierr; 10c4762a1bSJed Brown PetscMPIInt rank,size; 11c4762a1bSJed Brown PetscInt N0=3,N1=3,N2=3,N=N0*N1*N2; 12c4762a1bSJed Brown PetscRandom rdm; 13c4762a1bSJed Brown PetscScalar a; 14c4762a1bSJed Brown PetscReal enorm; 15c4762a1bSJed Brown Vec x,y,z,input,output; 16c4762a1bSJed Brown PetscBool view=PETSC_FALSE,use_interface=PETSC_TRUE; 17c4762a1bSJed Brown Mat A; 18c4762a1bSJed Brown PetscInt DIM, dim[3],vsize; 19c4762a1bSJed Brown PetscReal fac; 20c4762a1bSJed Brown 21c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 22c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 23c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers"); 24c4762a1bSJed Brown #endif 25c4762a1bSJed Brown 26ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRMPI(ierr); 27ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRMPI(ierr); 28c4762a1bSJed Brown 29c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD, &rdm);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr); 31c4762a1bSJed Brown 32c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&input);CHKERRQ(ierr); 33c4762a1bSJed Brown ierr = VecSetSizes(input,PETSC_DECIDE,N);CHKERRQ(ierr); 34c4762a1bSJed Brown ierr = VecSetFromOptions(input);CHKERRQ(ierr); 35c4762a1bSJed Brown ierr = VecSetRandom(input,rdm);CHKERRQ(ierr); 361e1ea65dSPierre Jolivet ierr = VecDuplicate(input,&output);CHKERRQ(ierr); 37c4762a1bSJed Brown /* ierr = VecGetSize(input,&vsize);CHKERRQ(ierr); */ 38c4762a1bSJed Brown /* printf("Size of the input Vector is %d\n",vsize); */ 39c4762a1bSJed Brown 40c4762a1bSJed Brown DIM = 3; 41c4762a1bSJed Brown dim[0] = N0; dim[1] = N1; dim[2] = N2; 42c4762a1bSJed Brown 43c4762a1bSJed Brown ierr = MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A);CHKERRQ(ierr); 44c4762a1bSJed Brown ierr = MatCreateVecs(A,&x,&y);CHKERRQ(ierr); 45c4762a1bSJed Brown ierr = MatCreateVecs(A,&z,NULL);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = VecGetSize(y,&vsize);CHKERRQ(ierr); 47c4762a1bSJed Brown printf("The vector size from the main routine is %d\n",vsize); 48c4762a1bSJed Brown 49c4762a1bSJed Brown ierr = InputTransformFFT(A,input,x);CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = MatMult(A,x,y);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = MatMultTranspose(A,y,z);CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = OutputTransformFFT(A,z,output);CHKERRQ(ierr); 53c4762a1bSJed Brown 54c4762a1bSJed Brown fac = 1.0/(PetscReal)N; 55c4762a1bSJed Brown ierr = VecScale(output,fac);CHKERRQ(ierr); 56c4762a1bSJed Brown 57c4762a1bSJed Brown ierr = VecAssemblyBegin(input);CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = VecAssemblyEnd(input);CHKERRQ(ierr); 59c4762a1bSJed Brown ierr = VecAssemblyBegin(output);CHKERRQ(ierr); 60c4762a1bSJed Brown ierr = VecAssemblyEnd(output);CHKERRQ(ierr); 61c4762a1bSJed Brown 62c4762a1bSJed Brown ierr = VecView(input,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 63c4762a1bSJed Brown ierr = VecView(output,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 64c4762a1bSJed Brown 65c4762a1bSJed Brown ierr = VecAXPY(output,-1.0,input);CHKERRQ(ierr); 66c4762a1bSJed Brown ierr = VecNorm(output,NORM_1,&enorm);CHKERRQ(ierr); 67c4762a1bSJed Brown /* if (enorm > 1.e-14) { */ 68*dd400576SPatrick Sanan if (rank == 0) { 69c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," Error norm of |x - z| %e\n",enorm);CHKERRQ(ierr); 70c4762a1bSJed Brown } 71c4762a1bSJed Brown /* } */ 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* ierr = MatCreateVecs(A,&z,NULL);CHKERRQ(ierr); */ 74c4762a1bSJed Brown /* printf("Vector size from ex148 %d\n",vsize); */ 75c4762a1bSJed Brown /* ierr = PetscObjectSetName((PetscObject) x, "Real space vector");CHKERRQ(ierr); */ 76c4762a1bSJed Brown /* ierr = PetscObjectSetName((PetscObject) y, "Frequency space vector");CHKERRQ(ierr); */ 77c4762a1bSJed Brown /* ierr = PetscObjectSetName((PetscObject) z, "Reconstructed vector");CHKERRQ(ierr); */ 78c4762a1bSJed Brown 79c4762a1bSJed Brown ierr = PetscFinalize(); 80c4762a1bSJed Brown return ierr; 81c4762a1bSJed Brown 82c4762a1bSJed Brown } 83c4762a1bSJed Brown 84