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