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