Actual source code: ex155.c
petsc-3.7.3 2016-08-01
1: static char help[]="This program illustrates the use of PETSc-fftw interface for parallel real DFT\n";
2: #include <petscmat.h>
3: #include <fftw3-mpi.h>
4: /*extern PetscErrorCode MatCreateVecsFFT(Mat,Vec *,Vec *,Vec *);*/
7: PetscInt main(PetscInt argc,char **args)
8: {
10: PetscMPIInt rank,size;
11: PetscInt N0=4096,N1=4096,N2=256,N3=10,N4=10,N=N0*N1;
12: PetscRandom rdm;
13: PetscReal enorm;
14: Vec x,y,z,input,output;
15: Mat A;
16: PetscInt DIM, dim[5],vsize,row,col;
17: PetscReal fac;
19: PetscInitialize(&argc,&args,(char*)0,help);
20: MPI_Comm_size(PETSC_COMM_WORLD, &size);
21: MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
23: #if defined(PETSC_USE_COMPLEX)
24: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "Example for Real DFT. Your current data type is complex!");
25: #endif
26: PetscRandomCreate(PETSC_COMM_WORLD, &rdm);
27: PetscRandomSetFromOptions(rdm);
28: VecCreate(PETSC_COMM_WORLD,&input);
29: VecSetSizes(input,PETSC_DECIDE,N);
30: VecSetFromOptions(input);
31: VecSetRandom(input,rdm);
32: VecDuplicate(input,&output);
34: DIM = 2; dim[0] = N0; dim[1] = N1; dim[2] = N2; dim[3] = N3; dim[4] = N4;
35: MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A);
36: MatGetLocalSize(A,&row,&col);
37: printf("The Matrix size is %d and %d from process %d\n",row,col,rank);
38: MatCreateVecsFFTW(A,&x,&y,&z);
40: VecGetSize(x,&vsize);
42: VecGetSize(z,&vsize);
43: printf("The vector size of output from the main routine is %d\n",vsize);
45: VecScatterPetscToFFTW(A,input,x);
46: /*VecDestroy(&input);*/
47: MatMult(A,x,y);
48: MatMultTranspose(A,y,z);
49: VecScatterFFTWToPetsc(A,z,output);
50: /*VecDestroy(&z);*/
51: fac = 1.0/(PetscReal)N;
52: VecScale(output,fac);
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);*/
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(&x);
69: VecDestroy(&y);
70: VecDestroy(&z);
71: VecDestroy(&output);
72: VecDestroy(&input);
73: MatDestroy(&A);
74: PetscRandomDestroy(&rdm);
75: PetscFinalize();
76: return 0;
78: }