xref: /petsc/src/mat/tests/ex157.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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*327415f7SBarry Smith   PetscFunctionBeginUser;
179566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
18c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
19c4762a1bSJed Brown   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers");
20c4762a1bSJed Brown #endif
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
23c4762a1bSJed Brown 
249566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
259566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
269566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&input));
279566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(input,PETSC_DECIDE,N));
289566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(input));
299566063dSJacob Faibussowitsch /*  PetscCall(VecSet(input,one)); */
309566063dSJacob Faibussowitsch /*  PetscCall(VecSetValue(input,1,two,INSERT_VALUES)); */
319566063dSJacob Faibussowitsch /*  PetscCall(VecSetValue(input,2,three,INSERT_VALUES)); */
329566063dSJacob Faibussowitsch /*  PetscCall(VecSetValue(input,3,three,INSERT_VALUES)); */
339566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(input,rdm));
349566063dSJacob Faibussowitsch /*  PetscCall(VecSetRandom(input,rdm)); */
359566063dSJacob Faibussowitsch /*  PetscCall(VecSetRandom(input,rdm)); */
369566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(input,&output));
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   DIM  = 2; dim[0] = N0; dim[1] = N1; dim[2] = N2; dim[3] = N3; dim[4] = N4;
399566063dSJacob Faibussowitsch   PetscCall(MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A));
409566063dSJacob Faibussowitsch   PetscCall(MatCreateVecsFFTW(A,&x,&y,&z));
419566063dSJacob Faibussowitsch /*  PetscCall(MatCreateVecs(A,&x,&y)); */
429566063dSJacob Faibussowitsch /*  PetscCall(MatCreateVecs(A,&z,NULL)); */
43c4762a1bSJed Brown 
449566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x,&vsize));
45c4762a1bSJed Brown   printf("The vector size  of input from the main routine is %d\n",vsize);
46c4762a1bSJed Brown 
479566063dSJacob Faibussowitsch   PetscCall(VecGetSize(z,&vsize));
48c4762a1bSJed Brown   printf("The vector size of output from the main routine is %d\n",vsize);
49c4762a1bSJed Brown 
509566063dSJacob Faibussowitsch   PetscCall(InputTransformFFT(A,input,x));
51c4762a1bSJed Brown 
529566063dSJacob Faibussowitsch   PetscCall(MatMult(A,x,y));
539566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(y));
549566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(y));
559566063dSJacob Faibussowitsch   PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
56c4762a1bSJed Brown 
579566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A,y,z));
58c4762a1bSJed Brown 
599566063dSJacob Faibussowitsch   PetscCall(OutputTransformFFT(A,z,output));
60c4762a1bSJed Brown   fac  = 1.0/(PetscReal)N;
619566063dSJacob Faibussowitsch   PetscCall(VecScale(output,fac));
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(input));
649566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(input));
659566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(output));
669566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(output));
67c4762a1bSJed Brown 
689566063dSJacob Faibussowitsch /*  PetscCall(VecView(input,PETSC_VIEWER_STDOUT_WORLD)); */
699566063dSJacob Faibussowitsch /*  PetscCall(VecView(output,PETSC_VIEWER_STDOUT_WORLD)); */
70c4762a1bSJed Brown 
719566063dSJacob Faibussowitsch   PetscCall(VecAXPY(output,-1.0,input));
729566063dSJacob Faibussowitsch   PetscCall(VecNorm(output,NORM_1,&enorm));
73c4762a1bSJed Brown /*  if (enorm > 1.e-14) { */
749566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"  Error norm of |x - z| %e\n",enorm));
75c4762a1bSJed Brown /*      } */
76c4762a1bSJed Brown 
779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&output));
789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&input));
799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&z));
829566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
839566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
849566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
85b122ec5aSJacob Faibussowitsch   return 0;
86c4762a1bSJed Brown 
87c4762a1bSJed Brown }
88