xref: /petsc/src/mat/tests/ex157.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[]="This program illustrates the use of PETSc-fftw interface for parallel real DFT\n";
2c4762a1bSJed Brown #include <petscmat.h>
3c4762a1bSJed Brown #include <fftw3-mpi.h>
4ad0e0ea3SStefano Zampini int main(int argc,char **args)
5c4762a1bSJed Brown {
6c4762a1bSJed Brown   PetscMPIInt    rank,size;
7c4762a1bSJed Brown   PetscInt       N0=2048,N1=2048,N2=3,N3=5,N4=5,N=N0*N1;
8c4762a1bSJed Brown   PetscRandom    rdm;
9c4762a1bSJed Brown   PetscReal      enorm;
10c4762a1bSJed Brown   Vec            x,y,z,input,output;
11c4762a1bSJed Brown   Mat            A;
12c4762a1bSJed Brown   PetscInt       DIM, dim[5],vsize;
13c4762a1bSJed Brown   PetscReal      fac;
14c4762a1bSJed Brown   PetscScalar    one=1,two=2,three=3;
15c4762a1bSJed Brown 
16*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
17c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
18c4762a1bSJed Brown   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers");
19c4762a1bSJed Brown #endif
205f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
215f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
22c4762a1bSJed Brown 
235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rdm));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&input));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(input,PETSC_DECIDE,N));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(input));
285f80ce2aSJacob Faibussowitsch /*  CHKERRQ(VecSet(input,one)); */
295f80ce2aSJacob Faibussowitsch /*  CHKERRQ(VecSetValue(input,1,two,INSERT_VALUES)); */
305f80ce2aSJacob Faibussowitsch /*  CHKERRQ(VecSetValue(input,2,three,INSERT_VALUES)); */
315f80ce2aSJacob Faibussowitsch /*  CHKERRQ(VecSetValue(input,3,three,INSERT_VALUES)); */
325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(input,rdm));
335f80ce2aSJacob Faibussowitsch /*  CHKERRQ(VecSetRandom(input,rdm)); */
345f80ce2aSJacob Faibussowitsch /*  CHKERRQ(VecSetRandom(input,rdm)); */
355f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(input,&output));
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   DIM  = 2; dim[0] = N0; dim[1] = N1; dim[2] = N2; dim[3] = N3; dim[4] = N4;
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecsFFTW(A,&x,&y,&z));
405f80ce2aSJacob Faibussowitsch /*  CHKERRQ(MatCreateVecs(A,&x,&y)); */
415f80ce2aSJacob Faibussowitsch /*  CHKERRQ(MatCreateVecs(A,&z,NULL)); */
42c4762a1bSJed Brown 
435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(x,&vsize));
44c4762a1bSJed Brown   printf("The vector size  of input from the main routine is %d\n",vsize);
45c4762a1bSJed Brown 
465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(z,&vsize));
47c4762a1bSJed Brown   printf("The vector size of output from the main routine is %d\n",vsize);
48c4762a1bSJed Brown 
495f80ce2aSJacob Faibussowitsch   CHKERRQ(InputTransformFFT(A,input,x));
50c4762a1bSJed Brown 
515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A,x,y));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(y));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(y));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
55c4762a1bSJed Brown 
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(A,y,z));
57c4762a1bSJed Brown 
585f80ce2aSJacob Faibussowitsch   CHKERRQ(OutputTransformFFT(A,z,output));
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) { */
735f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %e\n",enorm));
74c4762a1bSJed Brown /*      } */
75c4762a1bSJed Brown 
765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&output));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&input));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&z));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rdm));
83*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
84*b122ec5aSJacob Faibussowitsch   return 0;
85c4762a1bSJed Brown 
86c4762a1bSJed Brown }
87