xref: /petsc/src/mat/tests/ex150.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown static char help[]="This program illustrates the use of PETSc-fftw interface for real DFT\n";
2c4762a1bSJed Brown #include <petscmat.h>
3c4762a1bSJed Brown #include <fftw3-mpi.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown extern PetscErrorCode InputTransformFFT(Mat,Vec,Vec);
6c4762a1bSJed Brown extern PetscErrorCode OutputTransformFFT(Mat,Vec,Vec);
7ad0e0ea3SStefano Zampini int main(int argc,char **args)
8c4762a1bSJed Brown {
9c4762a1bSJed Brown   PetscErrorCode ierr;
10c4762a1bSJed Brown   PetscMPIInt    rank,size;
11c4762a1bSJed Brown   PetscInt       N0=3,N1=3,N2=3,N3=3,N4=3,N=N0*N1*N2*N3;
12c4762a1bSJed Brown   PetscRandom    rdm;
13c4762a1bSJed Brown   PetscReal      enorm;
14c4762a1bSJed Brown   Vec            x,y,z,input,output;
15c4762a1bSJed Brown   Mat            A;
16c4762a1bSJed Brown   PetscInt       DIM, dim[5],vsize;
17c4762a1bSJed Brown   PetscReal      fac;
18c4762a1bSJed Brown   PetscScalar    one=1;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
21c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
22c4762a1bSJed Brown   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers");
23c4762a1bSJed Brown #endif
24c4762a1bSJed Brown 
25*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
26*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
27c4762a1bSJed Brown 
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rdm));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&input));
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(input,PETSC_DECIDE,N));
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(input));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(input,rdm));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(input));
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(input));
36*5f80ce2aSJacob Faibussowitsch /*  CHKERRQ(VecView(input,PETSC_VIEWER_STDOUT_WORLD)); */
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(input,&output));
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   DIM    = 4;
40c4762a1bSJed Brown   dim[0] = N0; dim[1] = N1; dim[2] = N2; dim[3] = N3; dim[4] = N4;
41c4762a1bSJed Brown 
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A));
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&x,&y));
44*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&z,NULL));
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(x,&vsize));
46c4762a1bSJed Brown   printf("The vector size from the main routine is %d\n",vsize);
47c4762a1bSJed Brown 
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(InputTransformFFT(A,input,x));
49c4762a1bSJed Brown 
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A,x,y));
51*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(y));
52*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(y));
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
54*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(A,y,z));
55*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(z,&vsize));
56c4762a1bSJed Brown   printf("The vector size of zfrom the main routine is %d\n",vsize);
57c4762a1bSJed Brown 
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(OutputTransformFFT(A,z,output));
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   fac  = 1.0/(PetscReal)N;
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(output,fac));
62c4762a1bSJed Brown 
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(input));
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(input));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(output));
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(output));
67c4762a1bSJed Brown 
68*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(input,PETSC_VIEWER_STDOUT_WORLD));
69*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(output,PETSC_VIEWER_STDOUT_WORLD));
70c4762a1bSJed Brown 
71*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(output,-1.0,input));
72*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(output,NORM_1,&enorm));
73c4762a1bSJed Brown /*  if (enorm > 1.e-14) { */
74dd400576SPatrick Sanan   if (rank == 0) {
75*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %e\n",enorm));
76c4762a1bSJed Brown   }
77c4762a1bSJed Brown /*      } */
78c4762a1bSJed Brown 
79*5f80ce2aSJacob Faibussowitsch /* CHKERRQ(MatCreateVecs(A,&z,NULL)); */
80c4762a1bSJed Brown /*  printf("Vector size from ex148 %d\n",vsize); */
81*5f80ce2aSJacob Faibussowitsch /*  CHKERRQ(PetscObjectSetName((PetscObject) x, "Real space vector")); */
82*5f80ce2aSJacob Faibussowitsch /*      CHKERRQ(PetscObjectSetName((PetscObject) y, "Frequency space vector")); */
83*5f80ce2aSJacob Faibussowitsch /*      CHKERRQ(PetscObjectSetName((PetscObject) z, "Reconstructed vector")); */
84c4762a1bSJed Brown 
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&output));
86*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&input));
87*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
89*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&z));
90*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
91*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rdm));
92c4762a1bSJed Brown   ierr = PetscFinalize();
93c4762a1bSJed Brown   return ierr;
94c4762a1bSJed Brown }
95