1c4762a1bSJed Brown static char help[] = "This program illustrates the use of PETSc-fftw interface for sequential real DFT\n";
2c4762a1bSJed Brown #include <petscmat.h>
3c4762a1bSJed Brown #include <fftw3-mpi.h>
main(int argc,char ** args)4d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
5d71ae5a4SJacob Faibussowitsch {
6c4762a1bSJed Brown PetscMPIInt rank, size;
7c4762a1bSJed Brown PetscInt N0 = 10, N1 = 10, N2 = 10, N3 = 10, N4 = 10, N = N0 * N1 * N2 * N3 * N4;
8c4762a1bSJed Brown PetscRandom rdm;
9c4762a1bSJed Brown PetscReal enorm;
10c4762a1bSJed Brown Vec x, y, z, input, output;
11c4762a1bSJed Brown Mat A;
12c4762a1bSJed Brown PetscInt DIM, dim[5], vsize;
13c4762a1bSJed Brown PetscReal fac;
14c4762a1bSJed Brown
15327415f7SBarry Smith PetscFunctionBeginUser;
16*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
17c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
18c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "This example requires real numbers");
19c4762a1bSJed Brown #endif
209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
22c4762a1bSJed Brown
23be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uni-processor example only");
249566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
259566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm));
269566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &input));
279566063dSJacob Faibussowitsch PetscCall(VecSetSizes(input, N, N));
289566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(input));
299566063dSJacob Faibussowitsch PetscCall(VecSetRandom(input, rdm));
309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(input, &output));
31c4762a1bSJed Brown
329371c9d4SSatish Balay DIM = 5;
339371c9d4SSatish Balay dim[0] = N0;
349371c9d4SSatish Balay dim[1] = N1;
359371c9d4SSatish Balay dim[2] = N2;
369371c9d4SSatish Balay dim[3] = N3;
379371c9d4SSatish Balay dim[4] = N4;
389566063dSJacob Faibussowitsch PetscCall(MatCreateFFT(PETSC_COMM_SELF, DIM, dim, MATFFTW, &A));
399566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, &y));
409566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &z, NULL));
41c4762a1bSJed Brown
429566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &vsize));
43c4762a1bSJed Brown printf("The vector size of input from the main routine is %d\n", vsize);
44c4762a1bSJed Brown
459566063dSJacob Faibussowitsch PetscCall(VecGetSize(z, &vsize));
46c4762a1bSJed Brown printf("The vector size of output from the main routine is %d\n", vsize);
47c4762a1bSJed Brown
489566063dSJacob Faibussowitsch PetscCall(InputTransformFFT(A, input, x));
49c4762a1bSJed Brown
509566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, y));
519566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, y, z));
52c4762a1bSJed Brown
539566063dSJacob Faibussowitsch PetscCall(OutputTransformFFT(A, z, output));
54c4762a1bSJed Brown fac = 1.0 / (PetscReal)N;
559566063dSJacob Faibussowitsch PetscCall(VecScale(output, fac));
56c4762a1bSJed Brown /*
579566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(input));
589566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(input));
599566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(output));
609566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(output));
61c4762a1bSJed Brown
629566063dSJacob Faibussowitsch PetscCall(VecView(input,PETSC_VIEWER_STDOUT_WORLD));
639566063dSJacob Faibussowitsch PetscCall(VecView(output,PETSC_VIEWER_STDOUT_WORLD));
64c4762a1bSJed Brown */
659566063dSJacob Faibussowitsch PetscCall(VecAXPY(output, -1.0, input));
669566063dSJacob Faibussowitsch PetscCall(VecNorm(output, NORM_1, &enorm));
67c4762a1bSJed Brown /* if (enorm > 1.e-14) { */
689566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Error norm of |x - z| %e\n", enorm));
69c4762a1bSJed Brown /* } */
70c4762a1bSJed Brown
719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&output));
729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&input));
739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y));
759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&z));
769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
779566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm));
789566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
79b122ec5aSJacob Faibussowitsch return 0;
80c4762a1bSJed Brown }
81