xref: /petsc/src/mat/tests/ex157.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode ierr;
7c4762a1bSJed Brown   PetscMPIInt    rank,size;
8c4762a1bSJed Brown   PetscInt       N0=2048,N1=2048,N2=3,N3=5,N4=5,N=N0*N1;
9c4762a1bSJed Brown   PetscRandom    rdm;
10c4762a1bSJed Brown   PetscReal      enorm;
11c4762a1bSJed Brown   Vec            x,y,z,input,output;
12c4762a1bSJed Brown   Mat            A;
13c4762a1bSJed Brown   PetscInt       DIM, dim[5],vsize;
14c4762a1bSJed Brown   PetscReal      fac;
15c4762a1bSJed Brown   PetscScalar    one=1,two=2,three=3;
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
18c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
19c4762a1bSJed Brown   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers");
20c4762a1bSJed Brown #endif
21*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
22*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
23c4762a1bSJed Brown 
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rdm));
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&input));
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(input,PETSC_DECIDE,N));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(input));
29*5f80ce2aSJacob Faibussowitsch /*  CHKERRQ(VecSet(input,one)); */
30*5f80ce2aSJacob Faibussowitsch /*  CHKERRQ(VecSetValue(input,1,two,INSERT_VALUES)); */
31*5f80ce2aSJacob Faibussowitsch /*  CHKERRQ(VecSetValue(input,2,three,INSERT_VALUES)); */
32*5f80ce2aSJacob Faibussowitsch /*  CHKERRQ(VecSetValue(input,3,three,INSERT_VALUES)); */
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(input,rdm));
34*5f80ce2aSJacob Faibussowitsch /*  CHKERRQ(VecSetRandom(input,rdm)); */
35*5f80ce2aSJacob Faibussowitsch /*  CHKERRQ(VecSetRandom(input,rdm)); */
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(input,&output));
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   DIM  = 2; dim[0] = N0; dim[1] = N1; dim[2] = N2; dim[3] = N3; dim[4] = N4;
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A));
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecsFFTW(A,&x,&y,&z));
41*5f80ce2aSJacob Faibussowitsch /*  CHKERRQ(MatCreateVecs(A,&x,&y)); */
42*5f80ce2aSJacob Faibussowitsch /*  CHKERRQ(MatCreateVecs(A,&z,NULL)); */
43c4762a1bSJed Brown 
44*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(x,&vsize));
45c4762a1bSJed Brown   printf("The vector size  of input from the main routine is %d\n",vsize);
46c4762a1bSJed Brown 
47*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(z,&vsize));
48c4762a1bSJed Brown   printf("The vector size of output from the main routine is %d\n",vsize);
49c4762a1bSJed Brown 
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(InputTransformFFT(A,input,x));
51c4762a1bSJed Brown 
52*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A,x,y));
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(y));
54*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(y));
55*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
56c4762a1bSJed Brown 
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(A,y,z));
58c4762a1bSJed Brown 
59*5f80ce2aSJacob Faibussowitsch   CHKERRQ(OutputTransformFFT(A,z,output));
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) { */
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %e\n",enorm));
75c4762a1bSJed Brown /*      } */
76c4762a1bSJed Brown 
77*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&output));
78*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&input));
79*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
80*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&z));
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
83*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rdm));
84c4762a1bSJed Brown   ierr = PetscFinalize();
85c4762a1bSJed Brown   return ierr;
86c4762a1bSJed Brown 
87c4762a1bSJed Brown }
88