xref: /petsc/src/mat/tests/ex148.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode ierr;
9c4762a1bSJed Brown   PetscMPIInt    rank,size;
10c4762a1bSJed Brown   PetscInt       N0=50,N1=50,N=N0*N1;
11c4762a1bSJed Brown   PetscRandom    rdm;
12c4762a1bSJed Brown   PetscReal      enorm;
13c4762a1bSJed Brown   Vec            x,y,z,input,output;
14c4762a1bSJed Brown   Mat            A;
15c4762a1bSJed Brown   PetscInt       DIM,dim[2];
16c4762a1bSJed Brown   PetscReal      fac;
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
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 
23*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
24*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   /* Create and set PETSc vectors 'input' and 'output' */
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rdm));
29c4762a1bSJed Brown 
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&input));
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(input,PETSC_DECIDE,N0*N1));
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(input));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(input,rdm));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(input,&output));
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)input, "Real space vector"));
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecsFFTW(A,&x,&y,&z));
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* Scatter PETSc vector 'input' to FFTW vector 'x' */
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterPetscToFFTW(A,input,x));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* Apply forward FFT */
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A,x,y));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   /* Apply backward FFT */
51*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(A,y,z));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* Scatter FFTW vector 'z' to PETSc vector 'output' */
54*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterFFTWToPetsc(A,z,output));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* Check accuracy */
57c4762a1bSJed Brown   fac  = 1.0/(PetscReal)N;
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(output,fac));
59*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(output,-1.0,input));
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(output,NORM_1,&enorm));
61dd400576SPatrick Sanan   if (enorm > 1.e-11 && rank == 0) {
62*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %e\n",enorm));
63c4762a1bSJed Brown   }
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* Free spaces */
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rdm));
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&input));
68*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&output));
69*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
70*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
71*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&z));
72*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   ierr = PetscFinalize();
75c4762a1bSJed Brown   return ierr;
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