xref: /petsc/src/mat/tests/ex150.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,N3=3,N4=3,N=N0*N1*N2*N3;
11c4762a1bSJed Brown   PetscRandom    rdm;
12c4762a1bSJed Brown   PetscReal      enorm;
13c4762a1bSJed Brown   Vec            x,y,z,input,output;
14c4762a1bSJed Brown   Mat            A;
15c4762a1bSJed Brown   PetscInt       DIM, dim[5],vsize;
16c4762a1bSJed Brown   PetscReal      fac;
17c4762a1bSJed Brown   PetscScalar    one=1;
18c4762a1bSJed Brown 
19*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
20c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
21c4762a1bSJed Brown   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers");
22c4762a1bSJed Brown #endif
23c4762a1bSJed Brown 
245f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
255f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
26c4762a1bSJed Brown 
275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rdm));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&input));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(input,PETSC_DECIDE,N));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(input));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(input,rdm));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(input));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(input));
355f80ce2aSJacob Faibussowitsch /*  CHKERRQ(VecView(input,PETSC_VIEWER_STDOUT_WORLD)); */
365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(input,&output));
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   DIM    = 4;
39c4762a1bSJed Brown   dim[0] = N0; dim[1] = N1; dim[2] = N2; dim[3] = N3; dim[4] = N4;
40c4762a1bSJed Brown 
415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&x,&y));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&z,NULL));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(x,&vsize));
45c4762a1bSJed Brown   printf("The vector size from the main routine is %d\n",vsize);
46c4762a1bSJed Brown 
475f80ce2aSJacob Faibussowitsch   CHKERRQ(InputTransformFFT(A,input,x));
48c4762a1bSJed Brown 
495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A,x,y));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(y));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(y));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(A,y,z));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(z,&vsize));
55c4762a1bSJed Brown   printf("The vector size of zfrom the main routine is %d\n",vsize);
56c4762a1bSJed Brown 
575f80ce2aSJacob Faibussowitsch   CHKERRQ(OutputTransformFFT(A,z,output));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   fac  = 1.0/(PetscReal)N;
605f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(output,fac));
61c4762a1bSJed Brown 
625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(input));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(input));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(output));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(output));
66c4762a1bSJed Brown 
675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(input,PETSC_VIEWER_STDOUT_WORLD));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(output,PETSC_VIEWER_STDOUT_WORLD));
69c4762a1bSJed Brown 
705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(output,-1.0,input));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(output,NORM_1,&enorm));
72c4762a1bSJed Brown /*  if (enorm > 1.e-14) { */
73dd400576SPatrick Sanan   if (rank == 0) {
745f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %e\n",enorm));
75c4762a1bSJed Brown   }
76c4762a1bSJed Brown /*      } */
77c4762a1bSJed Brown 
785f80ce2aSJacob Faibussowitsch /* CHKERRQ(MatCreateVecs(A,&z,NULL)); */
79c4762a1bSJed Brown /*  printf("Vector size from ex148 %d\n",vsize); */
805f80ce2aSJacob Faibussowitsch /*  CHKERRQ(PetscObjectSetName((PetscObject) x, "Real space vector")); */
815f80ce2aSJacob Faibussowitsch /*      CHKERRQ(PetscObjectSetName((PetscObject) y, "Frequency space vector")); */
825f80ce2aSJacob Faibussowitsch /*      CHKERRQ(PetscObjectSetName((PetscObject) z, "Reconstructed vector")); */
83c4762a1bSJed Brown 
845f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&output));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&input));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&z));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rdm));
91*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
92*b122ec5aSJacob Faibussowitsch   return 0;
93c4762a1bSJed Brown }
94