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> 4c4762a1bSJed Brown /*extern PetscErrorCode MatCreateVecsFFT(Mat,Vec *,Vec *,Vec *);*/ 59371c9d4SSatish Balay int main(int argc, char **args) { 6c4762a1bSJed Brown PetscMPIInt rank, size; 7c4762a1bSJed Brown PetscInt N0 = 4096, N1 = 4096, N2 = 256, N3 = 10, N4 = 10, 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, row, col; 13c4762a1bSJed Brown PetscReal fac; 14c4762a1bSJed Brown 15327415f7SBarry Smith PetscFunctionBeginUser; 169566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 19c4762a1bSJed Brown 20c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 21c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Example for Real DFT. Your current data type is complex!"); 22c4762a1bSJed Brown #endif 239566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm)); 249566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 259566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &input)); 269566063dSJacob Faibussowitsch PetscCall(VecSetSizes(input, PETSC_DECIDE, N)); 279566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(input)); 289566063dSJacob Faibussowitsch PetscCall(VecSetRandom(input, rdm)); 299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(input, &output)); 30c4762a1bSJed Brown 319371c9d4SSatish Balay DIM = 2; 329371c9d4SSatish Balay dim[0] = N0; 339371c9d4SSatish Balay dim[1] = N1; 349371c9d4SSatish Balay dim[2] = N2; 359371c9d4SSatish Balay dim[3] = N3; 369371c9d4SSatish Balay dim[4] = N4; 379566063dSJacob Faibussowitsch PetscCall(MatCreateFFT(PETSC_COMM_WORLD, DIM, dim, MATFFTW, &A)); 389566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &row, &col)); 39c4762a1bSJed Brown printf("The Matrix size is %d and %d from process %d\n", row, col, rank); 409566063dSJacob Faibussowitsch PetscCall(MatCreateVecsFFTW(A, &x, &y, &z)); 41c4762a1bSJed Brown 429566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &vsize)); 43c4762a1bSJed Brown 449566063dSJacob Faibussowitsch PetscCall(VecGetSize(z, &vsize)); 45c4762a1bSJed Brown printf("The vector size of output from the main routine is %d\n", vsize); 46c4762a1bSJed Brown 479566063dSJacob Faibussowitsch PetscCall(VecScatterPetscToFFTW(A, input, x)); 489566063dSJacob Faibussowitsch /*PetscCall(VecDestroy(&input));*/ 499566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, y)); 509566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, y, z)); 519566063dSJacob Faibussowitsch PetscCall(VecScatterFFTWToPetsc(A, z, output)); 529566063dSJacob Faibussowitsch /*PetscCall(VecDestroy(&z));*/ 53c4762a1bSJed Brown fac = 1.0 / (PetscReal)N; 549566063dSJacob Faibussowitsch PetscCall(VecScale(output, fac)); 55c4762a1bSJed Brown 569566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(input)); 579566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(input)); 589566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(output)); 599566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(output)); 60c4762a1bSJed Brown 619566063dSJacob Faibussowitsch /* PetscCall(VecView(input,PETSC_VIEWER_STDOUT_WORLD));*/ 629566063dSJacob Faibussowitsch /* PetscCall(VecView(output,PETSC_VIEWER_STDOUT_WORLD));*/ 63c4762a1bSJed Brown 649566063dSJacob Faibussowitsch PetscCall(VecAXPY(output, -1.0, input)); 659566063dSJacob Faibussowitsch PetscCall(VecNorm(output, NORM_1, &enorm)); 66*48a46eb9SPierre Jolivet if (enorm > 1.e-14) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Error norm of |x - z| %e\n", enorm)); 67c4762a1bSJed Brown 689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&z)); 719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&output)); 729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&input)); 739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 749566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 759566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 76b122ec5aSJacob Faibussowitsch return 0; 77c4762a1bSJed Brown } 78