Actual source code: ex148.c
petsc-3.8.4 2018-03-24
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>
5: int main(int argc,char **args)
6: {
8: PetscMPIInt rank,size;
9: PetscInt N0=50,N1=50,N=N0*N1;
10: PetscRandom rdm;
11: PetscReal enorm;
12: Vec x,y,z,input,output;
13: Mat A;
14: PetscInt DIM,dim[2];
15: PetscReal fac;
17: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
18: #if defined(PETSC_USE_COMPLEX)
19: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers");
20: #endif
22: MPI_Comm_size(PETSC_COMM_WORLD, &size);
23: MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
25: /* Create and set PETSc vectors 'input' and 'output' */
26: PetscRandomCreate(PETSC_COMM_WORLD, &rdm);
27: PetscRandomSetFromOptions(rdm);
29: VecCreate(PETSC_COMM_WORLD,&input);
30: VecSetSizes(input,PETSC_DECIDE,N0*N1);
31: VecSetFromOptions(input);
32: VecSetRandom(input,rdm);
33: VecDuplicate(input,&output);
34: PetscObjectSetName((PetscObject)input, "Real space vector");
35: PetscObjectSetName((PetscObject)output, "Reconstructed vector");
37: /* Get FFTW vectors 'x', 'y' and 'z' */
38: DIM = 2;
39: dim[0] = N0; dim[1] = N1;
40: MatCreateFFT(PETSC_COMM_WORLD,DIM,dim,MATFFTW,&A);
41: MatCreateVecsFFTW(A,&x,&y,&z);
43: /* Scatter PETSc vector 'input' to FFTW vector 'x' */
44: VecScatterPetscToFFTW(A,input,x);
46: /* Apply forward FFT */
47: MatMult(A,x,y);
49: /* Apply backward FFT */
50: MatMultTranspose(A,y,z);
52: /* Scatter FFTW vector 'z' to PETSc vector 'output' */
53: VecScatterFFTWToPetsc(A,z,output);
55: /* Check accuracy */
56: fac = 1.0/(PetscReal)N;
57: VecScale(output,fac);
58: VecAXPY(output,-1.0,input);
59: VecNorm(output,NORM_1,&enorm);
60: if (enorm > 1.e-11 && !rank) {
61: PetscPrintf(PETSC_COMM_SELF," Error norm of |x - z| %e\n",enorm);
62: }
64: /* Free spaces */
65: PetscRandomDestroy(&rdm);
66: VecDestroy(&input);
67: VecDestroy(&output);
68: VecDestroy(&x);
69: VecDestroy(&y);
70: VecDestroy(&z);
71: MatDestroy(&A);
73: PetscFinalize();
74: return ierr;
75: }