Actual source code: ex153.c
petsc-3.5.4 2015-05-23
1: static char help[]="This program illustrates the use of PETSc-fftw interface for sequential real DFT\n";
2: #include <petscmat.h>
3: #include <fftw3-mpi.h>
6: PetscInt main(PetscInt argc,char **args)
7: {
9: PetscMPIInt rank,size;
10: PetscInt N0=10,N1=10,N2=10,N3=10,N4=10,N=N0*N1*N2*N3*N4;
11: PetscRandom rdm;
12: PetscReal enorm;
13: Vec x,y,z,input,output;
14: Mat A;
15: PetscInt DIM, dim[5],vsize;
16: PetscReal fac;
18: PetscInitialize(&argc,&args,(char*)0,help);
19: #if defined(PETSC_USE_COMPLEX)
20: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers");
21: #endif
22: MPI_Comm_size(PETSC_COMM_WORLD, &size);
23: MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
25: if (size!=1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This is a uni-processor example only");
26: PetscRandomCreate(PETSC_COMM_SELF, &rdm);
27: PetscRandomSetFromOptions(rdm);
28: VecCreate(PETSC_COMM_SELF,&input);
29: VecSetSizes(input,N,N);
30: VecSetFromOptions(input);
31: VecSetRandom(input,rdm);
32: VecDuplicate(input,&output);
34: DIM = 5; dim[0] = N0; dim[1] = N1; dim[2] = N2; dim[3] = N3; dim[4] = N4;
35: MatCreateFFT(PETSC_COMM_SELF,DIM,dim,MATFFTW,&A);
36: MatGetVecs(A,&x,&y);
37: MatGetVecs(A,&z,NULL);
39: VecGetSize(x,&vsize);
40: printf("The vector size of input from the main routine is %d\n",vsize);
42: VecGetSize(z,&vsize);
43: printf("The vector size of output from the main routine is %d\n",vsize);
45: InputTransformFFT(A,input,x);
47: MatMult(A,x,y);
48: MatMultTranspose(A,y,z);
50: OutputTransformFFT(A,z,output);
51: fac = 1.0/(PetscReal)N;
52: VecScale(output,fac);
53: /*
54: VecAssemblyBegin(input);
55: VecAssemblyEnd(input);
56: VecAssemblyBegin(output);
57: VecAssemblyEnd(output);
59: VecView(input,PETSC_VIEWER_STDOUT_WORLD);
60: VecView(output,PETSC_VIEWER_STDOUT_WORLD);
61: */
62: VecAXPY(output,-1.0,input);
63: VecNorm(output,NORM_1,&enorm);
64: /* if (enorm > 1.e-14) { */
65: PetscPrintf(PETSC_COMM_SELF," Error norm of |x - z| %e\n",enorm);
66: /* } */
68: VecDestroy(&output);
69: VecDestroy(&input);
70: VecDestroy(&x);
71: VecDestroy(&y);
72: VecDestroy(&z);
73: MatDestroy(&A);
74: PetscRandomDestroy(&rdm);
75: PetscFinalize();
76: return 0;
78: }