xref: /petsc/src/mat/tests/ex148.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
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 
main(int argc,char ** args)6d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   PetscMPIInt rank, size;
9c4762a1bSJed Brown   PetscInt    N0 = 50, N1 = 50, N = N0 * N1;
10c4762a1bSJed Brown   PetscRandom rdm;
11c4762a1bSJed Brown   PetscReal   enorm;
12c4762a1bSJed Brown   Vec         x, y, z, input, output;
13c4762a1bSJed Brown   Mat         A;
14c4762a1bSJed Brown   PetscInt    DIM, dim[2];
15c4762a1bSJed Brown   PetscReal   fac;
16c4762a1bSJed Brown 
17327415f7SBarry Smith   PetscFunctionBeginUser;
18c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
19c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
20c4762a1bSJed Brown   SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "This example requires real numbers");
21c4762a1bSJed Brown #endif
22c4762a1bSJed Brown 
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   /* Create and set PETSc vectors 'input' and 'output' */
279566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
289566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
29c4762a1bSJed Brown 
309566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &input));
319566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(input, PETSC_DECIDE, N0 * N1));
329566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(input));
339566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(input, rdm));
349566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(input, &output));
359566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)input, "Real space vector"));
369566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)output, "Reconstructed vector"));
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   /* Get FFTW vectors 'x', 'y' and 'z' */
39c4762a1bSJed Brown   DIM    = 2;
409371c9d4SSatish Balay   dim[0] = N0;
419371c9d4SSatish Balay   dim[1] = N1;
429566063dSJacob Faibussowitsch   PetscCall(MatCreateFFT(PETSC_COMM_WORLD, DIM, dim, MATFFTW, &A));
439566063dSJacob Faibussowitsch   PetscCall(MatCreateVecsFFTW(A, &x, &y, &z));
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* Scatter PETSc vector 'input' to FFTW vector 'x' */
469566063dSJacob Faibussowitsch   PetscCall(VecScatterPetscToFFTW(A, input, x));
47c4762a1bSJed Brown 
48c4762a1bSJed Brown   /* Apply forward FFT */
499566063dSJacob Faibussowitsch   PetscCall(MatMult(A, x, y));
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* Apply backward FFT */
529566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A, y, z));
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   /* Scatter FFTW vector 'z' to PETSc vector 'output' */
559566063dSJacob Faibussowitsch   PetscCall(VecScatterFFTWToPetsc(A, z, output));
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   /* Check accuracy */
58c4762a1bSJed Brown   fac = 1.0 / (PetscReal)N;
599566063dSJacob Faibussowitsch   PetscCall(VecScale(output, fac));
609566063dSJacob Faibussowitsch   PetscCall(VecAXPY(output, -1.0, input));
619566063dSJacob Faibussowitsch   PetscCall(VecNorm(output, NORM_1, &enorm));
6248a46eb9SPierre Jolivet   if (enorm > 1.e-11 && rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Error norm of |x - z| %e\n", enorm));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   /* Free spaces */
659566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&input));
679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&output));
689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&z));
719566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
72c4762a1bSJed Brown 
739566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
74b122ec5aSJacob Faibussowitsch   return 0;
75c4762a1bSJed Brown }
76c4762a1bSJed Brown 
77c4762a1bSJed Brown /*TEST
78c4762a1bSJed Brown 
79c4762a1bSJed Brown    build:
80c4762a1bSJed Brown       requires: fftw !complex
81c4762a1bSJed Brown 
82c4762a1bSJed Brown    test:
83*3886731fSPierre Jolivet       output_file: output/empty.out
84c4762a1bSJed Brown 
85c4762a1bSJed Brown    test:
86c4762a1bSJed Brown       suffix: 2
87c4762a1bSJed Brown       nsize: 3
88*3886731fSPierre Jolivet       output_file: output/empty.out
89c4762a1bSJed Brown 
90c4762a1bSJed Brown TEST*/
91