xref: /petsc/src/mat/tests/ex155.c (revision 6a98f8dc3f2c9149905a87dc2e9d0fedaf64e09a)
1 static char help[]="This program illustrates the use of PETSc-fftw interface for parallel real DFT\n";
2 #include <petscmat.h>
3 #include <fftw3-mpi.h>
4 /*extern PetscErrorCode MatCreateVecsFFT(Mat,Vec *,Vec *,Vec *);*/
5 PetscInt main(PetscInt argc,char **args)
6 {
7   PetscErrorCode ierr;
8   PetscMPIInt    rank,size;
9   PetscInt       N0=4096,N1=4096,N2=256,N3=10,N4=10,N=N0*N1;
10   PetscRandom    rdm;
11   PetscReal      enorm;
12   Vec            x,y,z,input,output;
13   Mat            A;
14   PetscInt       DIM, dim[5],vsize,row,col;
15   PetscReal      fac;
16 
17   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
18   ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRQ(ierr);
19   ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr);
20 
21 #if defined(PETSC_USE_COMPLEX)
22   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "Example for Real DFT. Your current data type is complex!");
23 #endif
24   ierr = PetscRandomCreate(PETSC_COMM_WORLD, &rdm);CHKERRQ(ierr);
25   ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr);
26   ierr = VecCreate(PETSC_COMM_WORLD,&input);CHKERRQ(ierr);
27   ierr = VecSetSizes(input,PETSC_DECIDE,N);CHKERRQ(ierr);
28   ierr = VecSetFromOptions(input);CHKERRQ(ierr);
29   ierr = VecSetRandom(input,rdm);CHKERRQ(ierr);
30   ierr = VecDuplicate(input,&output);
31 
32   DIM  = 2; dim[0] = N0; dim[1] = N1; dim[2] = N2; dim[3] = N3; dim[4] = N4;
33   ierr = MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A);CHKERRQ(ierr);
34   ierr = MatGetLocalSize(A,&row,&col);CHKERRQ(ierr);
35   printf("The Matrix size  is %d and %d from process %d\n",row,col,rank);
36   ierr = MatCreateVecsFFTW(A,&x,&y,&z);CHKERRQ(ierr);
37 
38   ierr = VecGetSize(x,&vsize);CHKERRQ(ierr);
39 
40   ierr = VecGetSize(z,&vsize);CHKERRQ(ierr);
41   printf("The vector size of output from the main routine is %d\n",vsize);
42 
43   ierr = VecScatterPetscToFFTW(A,input,x);CHKERRQ(ierr);
44   /*ierr = VecDestroy(&input);CHKERRQ(ierr);*/
45   ierr = MatMult(A,x,y);CHKERRQ(ierr);
46   ierr = MatMultTranspose(A,y,z);CHKERRQ(ierr);
47   ierr = VecScatterFFTWToPetsc(A,z,output);CHKERRQ(ierr);
48   /*ierr = VecDestroy(&z);CHKERRQ(ierr);*/
49   fac  = 1.0/(PetscReal)N;
50   ierr = VecScale(output,fac);CHKERRQ(ierr);
51 
52   ierr = VecAssemblyBegin(input);CHKERRQ(ierr);
53   ierr = VecAssemblyEnd(input);CHKERRQ(ierr);
54   ierr = VecAssemblyBegin(output);CHKERRQ(ierr);
55   ierr = VecAssemblyEnd(output);CHKERRQ(ierr);
56 
57 /*  ierr = VecView(input,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/
58 /*  ierr = VecView(output,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/
59 
60   ierr = VecAXPY(output,-1.0,input);CHKERRQ(ierr);
61   ierr = VecNorm(output,NORM_1,&enorm);CHKERRQ(ierr);
62   if (enorm > 1.e-14) {
63     ierr = PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %e\n",enorm);CHKERRQ(ierr);
64   }
65 
66   ierr = VecDestroy(&x);CHKERRQ(ierr);
67   ierr = VecDestroy(&y);CHKERRQ(ierr);
68   ierr = VecDestroy(&z);CHKERRQ(ierr);
69   ierr = VecDestroy(&output);CHKERRQ(ierr);
70   ierr = VecDestroy(&input);CHKERRQ(ierr);
71   ierr = MatDestroy(&A);CHKERRQ(ierr);
72   ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr);
73   ierr = PetscFinalize();
74   return ierr;
75 
76 }
77