xref: /petsc/src/mat/tests/ex148.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown static char help[] = "This program illustrates the use of PETSc-fftw interface for real 2D DFT.\n\
2c4762a1bSJed Brown                     See ~petsc/src/mat/tests/ex158.c for general cases. \n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
69371c9d4SSatish Balay int main(int argc, char **args) {
7c4762a1bSJed Brown   PetscMPIInt rank, size;
8c4762a1bSJed Brown   PetscInt    N0 = 50, N1 = 50, N = N0 * N1;
9c4762a1bSJed Brown   PetscRandom rdm;
10c4762a1bSJed Brown   PetscReal   enorm;
11c4762a1bSJed Brown   Vec         x, y, z, input, output;
12c4762a1bSJed Brown   Mat         A;
13c4762a1bSJed Brown   PetscInt    DIM, dim[2];
14c4762a1bSJed Brown   PetscReal   fac;
15c4762a1bSJed Brown 
16327415f7SBarry 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
21c4762a1bSJed Brown 
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
24c4762a1bSJed Brown 
25c4762a1bSJed Brown   /* Create and set PETSc vectors 'input' and 'output' */
269566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
279566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
28c4762a1bSJed Brown 
299566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &input));
309566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(input, PETSC_DECIDE, N0 * N1));
319566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(input));
329566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(input, rdm));
339566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(input, &output));
349566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)input, "Real space vector"));
359566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)output, "Reconstructed vector"));
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /* Get FFTW vectors 'x', 'y' and 'z' */
38c4762a1bSJed Brown   DIM    = 2;
399371c9d4SSatish Balay   dim[0] = N0;
409371c9d4SSatish Balay   dim[1] = N1;
419566063dSJacob Faibussowitsch   PetscCall(MatCreateFFT(PETSC_COMM_WORLD, DIM, dim, MATFFTW, &A));
429566063dSJacob Faibussowitsch   PetscCall(MatCreateVecsFFTW(A, &x, &y, &z));
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* Scatter PETSc vector 'input' to FFTW vector 'x' */
459566063dSJacob Faibussowitsch   PetscCall(VecScatterPetscToFFTW(A, input, x));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* Apply forward FFT */
489566063dSJacob Faibussowitsch   PetscCall(MatMult(A, x, y));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   /* Apply backward FFT */
519566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A, y, z));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* Scatter FFTW vector 'z' to PETSc vector 'output' */
549566063dSJacob Faibussowitsch   PetscCall(VecScatterFFTWToPetsc(A, z, output));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* Check accuracy */
57c4762a1bSJed Brown   fac = 1.0 / (PetscReal)N;
589566063dSJacob Faibussowitsch   PetscCall(VecScale(output, fac));
599566063dSJacob Faibussowitsch   PetscCall(VecAXPY(output, -1.0, input));
609566063dSJacob Faibussowitsch   PetscCall(VecNorm(output, NORM_1, &enorm));
61*48a46eb9SPierre Jolivet   if (enorm > 1.e-11 && rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Error norm of |x - z| %e\n", enorm));
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   /* Free spaces */
649566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&input));
669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&output));
679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&z));
709566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
73b122ec5aSJacob Faibussowitsch   return 0;
74c4762a1bSJed Brown }
75c4762a1bSJed Brown 
76c4762a1bSJed Brown /*TEST
77c4762a1bSJed Brown 
78c4762a1bSJed Brown    build:
79c4762a1bSJed Brown       requires: fftw !complex
80c4762a1bSJed Brown 
81c4762a1bSJed Brown    test:
82c4762a1bSJed Brown       output_file: output/ex148.out
83c4762a1bSJed Brown 
84c4762a1bSJed Brown    test:
85c4762a1bSJed Brown       suffix: 2
86c4762a1bSJed Brown       nsize: 3
87c4762a1bSJed Brown       output_file: output/ex148.out
88c4762a1bSJed Brown 
89c4762a1bSJed Brown TEST*/
90