Actual source code: ex148.c
petsc-3.7.3 2016-08-01
1: static char help[]="This program illustrates the use of PETSc-fftw interface for real 2D DFT.\n\
2: See ~petsc/src/mat/examples/tests/ex158.c for general cases. \n\n";
3: #include <petscmat.h>
7: int main(int argc,char **args)
8: {
10: PetscMPIInt rank,size;
11: PetscInt N0=50,N1=50,N=N0*N1;
12: PetscRandom rdm;
13: PetscReal enorm;
14: Vec x,y,z,input,output;
15: Mat A;
16: PetscInt DIM,dim[2];
17: PetscReal fac;
19: PetscInitialize(&argc,&args,(char*)0,help);
20: #if defined(PETSC_USE_COMPLEX)
21: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers");
22: #endif
24: MPI_Comm_size(PETSC_COMM_WORLD, &size);
25: MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
27: /* Create and set PETSc vectors 'input' and 'output' */
28: PetscRandomCreate(PETSC_COMM_WORLD, &rdm);
29: PetscRandomSetFromOptions(rdm);
31: VecCreate(PETSC_COMM_WORLD,&input);
32: VecSetSizes(input,PETSC_DECIDE,N0*N1);
33: VecSetFromOptions(input);
34: VecSetRandom(input,rdm);
35: VecDuplicate(input,&output);
36: PetscObjectSetName((PetscObject)input, "Real space vector");
37: PetscObjectSetName((PetscObject)output, "Reconstructed vector");
39: /* Get FFTW vectors 'x', 'y' and 'z' */
40: DIM = 2;
41: dim[0] = N0; dim[1] = N1;
42: MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A);
43: MatCreateVecsFFTW(A,&x,&y,&z);
45: /* Scatter PETSc vector 'input' to FFTW vector 'x' */
46: VecScatterPetscToFFTW(A,input,x);
48: /* Apply forward FFT */
49: MatMult(A,x,y);
51: /* Apply backward FFT */
52: MatMultTranspose(A,y,z);
54: /* Scatter FFTW vector 'z' to PETSc vector 'output' */
55: VecScatterFFTWToPetsc(A,z,output);
57: /* Check accuracy */
58: fac = 1.0/(PetscReal)N;
59: VecScale(output,fac);
60: VecAXPY(output,-1.0,input);
61: VecNorm(output,NORM_1,&enorm);
62: if (enorm > 1.e-11 && !rank) {
63: PetscPrintf(PETSC_COMM_SELF," Error norm of |x - z| %e\n",enorm);
64: }
66: /* Free spaces */
67: PetscRandomDestroy(&rdm);
68: VecDestroy(&input);
69: VecDestroy(&output);
70: VecDestroy(&x);
71: VecDestroy(&y);
72: VecDestroy(&z);
73: MatDestroy(&A);
75: PetscFinalize();
76: return 0;
77: }