xref: /petsc/src/mat/tests/ex149.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,N=N0*N1*N2;
12c4762a1bSJed Brown   PetscRandom    rdm;
13c4762a1bSJed Brown   PetscScalar    a;
14c4762a1bSJed Brown   PetscReal      enorm;
15c4762a1bSJed Brown   Vec            x,y,z,input,output;
16c4762a1bSJed Brown   PetscBool      view=PETSC_FALSE,use_interface=PETSC_TRUE;
17c4762a1bSJed Brown   Mat            A;
18c4762a1bSJed Brown   PetscInt       DIM, dim[3],vsize;
19c4762a1bSJed Brown   PetscReal      fac;
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
22c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
23c4762a1bSJed Brown   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers");
24c4762a1bSJed Brown #endif
25c4762a1bSJed Brown 
26*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
27*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
28c4762a1bSJed Brown 
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rdm));
31c4762a1bSJed Brown 
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&input));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(input,PETSC_DECIDE,N));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(input));
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(input,rdm));
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(input,&output));
37*5f80ce2aSJacob Faibussowitsch /*  CHKERRQ(VecGetSize(input,&vsize)); */
38c4762a1bSJed Brown /*  printf("Size of the input Vector is %d\n",vsize); */
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   DIM    = 3;
41c4762a1bSJed Brown   dim[0] = N0; dim[1] = N1; dim[2] = N2;
42c4762a1bSJed Brown 
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A));
44*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&x,&y));
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&z,NULL));
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(y,&vsize));
47c4762a1bSJed Brown   printf("The vector size from the main routine is %d\n",vsize);
48c4762a1bSJed Brown 
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(InputTransformFFT(A,input,x));
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A,x,y));
51*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(A,y,z));
52*5f80ce2aSJacob Faibussowitsch   CHKERRQ(OutputTransformFFT(A,z,output));
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   fac  = 1.0/(PetscReal)N;
55*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(output,fac));
56c4762a1bSJed Brown 
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(input));
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(input));
59*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(output));
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(output));
61c4762a1bSJed Brown 
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(input,PETSC_VIEWER_STDOUT_WORLD));
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(output,PETSC_VIEWER_STDOUT_WORLD));
64c4762a1bSJed Brown 
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(output,-1.0,input));
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(output,NORM_1,&enorm));
67c4762a1bSJed Brown /*  if (enorm > 1.e-14) { */
68dd400576SPatrick Sanan   if (rank == 0) {
69*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %e\n",enorm));
70c4762a1bSJed Brown   }
71c4762a1bSJed Brown /*      } */
72c4762a1bSJed Brown 
73*5f80ce2aSJacob Faibussowitsch /* CHKERRQ(MatCreateVecs(A,&z,NULL)); */
74c4762a1bSJed Brown /*  printf("Vector size from ex148 %d\n",vsize); */
75*5f80ce2aSJacob Faibussowitsch /*  CHKERRQ(PetscObjectSetName((PetscObject) x, "Real space vector")); */
76*5f80ce2aSJacob Faibussowitsch /*      CHKERRQ(PetscObjectSetName((PetscObject) y, "Frequency space vector")); */
77*5f80ce2aSJacob Faibussowitsch /*      CHKERRQ(PetscObjectSetName((PetscObject) z, "Reconstructed vector")); */
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   ierr = PetscFinalize();
80c4762a1bSJed Brown   return ierr;
81c4762a1bSJed Brown 
82c4762a1bSJed Brown }
83