xref: /petsc/src/mat/tests/ex148.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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*327415f7SBarry Smith   PetscFunctionBeginUser;
189566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
19c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
20c4762a1bSJed Brown   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers");
21c4762a1bSJed Brown #endif
22c4762a1bSJed Brown 
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   /* Create and set PETSc vectors 'input' and 'output' */
279566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
289566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
29c4762a1bSJed Brown 
309566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&input));
319566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(input,PETSC_DECIDE,N0*N1));
329566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(input));
339566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(input,rdm));
349566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(input,&output));
359566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)input, "Real space vector"));
369566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)output, "Reconstructed vector"));
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   /* Get FFTW vectors 'x', 'y' and 'z' */
39c4762a1bSJed Brown   DIM    = 2;
40c4762a1bSJed Brown   dim[0] = N0; dim[1] = N1;
419566063dSJacob Faibussowitsch   PetscCall(MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A));
429566063dSJacob Faibussowitsch   PetscCall(MatCreateVecsFFTW(A,&x,&y,&z));
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* Scatter PETSc vector 'input' to FFTW vector 'x' */
459566063dSJacob Faibussowitsch   PetscCall(VecScatterPetscToFFTW(A,input,x));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* Apply forward FFT */
489566063dSJacob Faibussowitsch   PetscCall(MatMult(A,x,y));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   /* Apply backward FFT */
519566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A,y,z));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* Scatter FFTW vector 'z' to PETSc vector 'output' */
549566063dSJacob Faibussowitsch   PetscCall(VecScatterFFTWToPetsc(A,z,output));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* Check accuracy */
57c4762a1bSJed Brown   fac  = 1.0/(PetscReal)N;
589566063dSJacob Faibussowitsch   PetscCall(VecScale(output,fac));
599566063dSJacob Faibussowitsch   PetscCall(VecAXPY(output,-1.0,input));
609566063dSJacob Faibussowitsch   PetscCall(VecNorm(output,NORM_1,&enorm));
61dd400576SPatrick Sanan   if (enorm > 1.e-11 && rank == 0) {
629566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %e\n",enorm));
63c4762a1bSJed Brown   }
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* Free spaces */
669566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&input));
689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&output));
699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&z));
729566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
73c4762a1bSJed Brown 
749566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
75b122ec5aSJacob Faibussowitsch   return 0;
76c4762a1bSJed Brown }
77c4762a1bSJed Brown 
78c4762a1bSJed Brown /*TEST
79c4762a1bSJed Brown 
80c4762a1bSJed Brown    build:
81c4762a1bSJed Brown       requires: fftw !complex
82c4762a1bSJed Brown 
83c4762a1bSJed Brown    test:
84c4762a1bSJed Brown       output_file: output/ex148.out
85c4762a1bSJed Brown 
86c4762a1bSJed Brown    test:
87c4762a1bSJed Brown       suffix: 2
88c4762a1bSJed Brown       nsize: 3
89c4762a1bSJed Brown       output_file: output/ex148.out
90c4762a1bSJed Brown 
91c4762a1bSJed Brown TEST*/
92