xref: /petsc/src/mat/tests/ex149.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscMPIInt    rank,size;
10c4762a1bSJed Brown   PetscInt       N0=3,N1=3,N2=3,N=N0*N1*N2;
11c4762a1bSJed Brown   PetscRandom    rdm;
12c4762a1bSJed Brown   PetscScalar    a;
13c4762a1bSJed Brown   PetscReal      enorm;
14c4762a1bSJed Brown   Vec            x,y,z,input,output;
15c4762a1bSJed Brown   PetscBool      view=PETSC_FALSE,use_interface=PETSC_TRUE;
16c4762a1bSJed Brown   Mat            A;
17c4762a1bSJed Brown   PetscInt       DIM, dim[3],vsize;
18c4762a1bSJed Brown   PetscReal      fac;
19c4762a1bSJed Brown 
20*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
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 
255f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
265f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
27c4762a1bSJed Brown 
285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rdm));
30c4762a1bSJed Brown 
315f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&input));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(input,PETSC_DECIDE,N));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(input));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(input,rdm));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(input,&output));
365f80ce2aSJacob Faibussowitsch /*  CHKERRQ(VecGetSize(input,&vsize)); */
37c4762a1bSJed Brown /*  printf("Size of the input Vector is %d\n",vsize); */
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   DIM    = 3;
40c4762a1bSJed Brown   dim[0] = N0; dim[1] = N1; dim[2] = N2;
41c4762a1bSJed Brown 
425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&x,&y));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&z,NULL));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(y,&vsize));
46c4762a1bSJed Brown   printf("The vector size from the main routine is %d\n",vsize);
47c4762a1bSJed Brown 
485f80ce2aSJacob Faibussowitsch   CHKERRQ(InputTransformFFT(A,input,x));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A,x,y));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(A,y,z));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(OutputTransformFFT(A,z,output));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   fac  = 1.0/(PetscReal)N;
545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(output,fac));
55c4762a1bSJed Brown 
565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(input));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(input));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(output));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(output));
60c4762a1bSJed Brown 
615f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(input,PETSC_VIEWER_STDOUT_WORLD));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(output,PETSC_VIEWER_STDOUT_WORLD));
63c4762a1bSJed Brown 
645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(output,-1.0,input));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(output,NORM_1,&enorm));
66c4762a1bSJed Brown /*  if (enorm > 1.e-14) { */
67dd400576SPatrick Sanan   if (rank == 0) {
685f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %e\n",enorm));
69c4762a1bSJed Brown   }
70c4762a1bSJed Brown /*      } */
71c4762a1bSJed Brown 
725f80ce2aSJacob Faibussowitsch /* CHKERRQ(MatCreateVecs(A,&z,NULL)); */
73c4762a1bSJed Brown /*  printf("Vector size from ex148 %d\n",vsize); */
745f80ce2aSJacob Faibussowitsch /*  CHKERRQ(PetscObjectSetName((PetscObject) x, "Real space vector")); */
755f80ce2aSJacob Faibussowitsch /*      CHKERRQ(PetscObjectSetName((PetscObject) y, "Frequency space vector")); */
765f80ce2aSJacob Faibussowitsch /*      CHKERRQ(PetscObjectSetName((PetscObject) z, "Reconstructed vector")); */
77c4762a1bSJed Brown 
78*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
79*b122ec5aSJacob Faibussowitsch   return 0;
80c4762a1bSJed Brown 
81c4762a1bSJed Brown }
82