1c4762a1bSJed Brown static char help[] = "Illustrate how to use mpi FFTW and PETSc-FFTW interface \n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown Usage:
5c4762a1bSJed Brown mpiexec -n <np> ./ex158 -use_FFTW_interface NO
6c4762a1bSJed Brown mpiexec -n <np> ./ex158 -use_FFTW_interface YES
7c4762a1bSJed Brown */
8c4762a1bSJed Brown
9c4762a1bSJed Brown #include <petscmat.h>
10c4762a1bSJed Brown #include <fftw3-mpi.h>
11c4762a1bSJed Brown
main(int argc,char ** args)12d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
13d71ae5a4SJacob Faibussowitsch {
14c4762a1bSJed Brown PetscMPIInt rank, size;
15c4762a1bSJed Brown PetscInt N0 = 50, N1 = 20, N = N0 * N1;
16c4762a1bSJed Brown PetscRandom rdm;
17c4762a1bSJed Brown PetscScalar a;
18c4762a1bSJed Brown PetscReal enorm;
19c4762a1bSJed Brown Vec x, y, z;
20c4762a1bSJed Brown PetscBool view = PETSC_FALSE, use_interface = PETSC_TRUE;
21c4762a1bSJed Brown
22327415f7SBarry Smith PetscFunctionBeginUser;
23*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
24c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
25c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "This example requires real numbers. Your current scalar type is complex");
26c4762a1bSJed Brown #endif
27c4762a1bSJed Brown
28d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex158");
299566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_FFTW_interface", "Use PETSc-FFTW interface", "ex158", use_interface, &use_interface, NULL));
30d0609cedSBarry Smith PetscOptionsEnd();
31c4762a1bSJed Brown
329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
34c4762a1bSJed Brown
359566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
369566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm));
37c4762a1bSJed Brown
38c4762a1bSJed Brown if (!use_interface) {
39c4762a1bSJed Brown /* Use mpi FFTW without PETSc-FFTW interface, 2D case only */
40c4762a1bSJed Brown /*---------------------------------------------------------*/
41c4762a1bSJed Brown fftw_plan fplan, bplan;
42c4762a1bSJed Brown fftw_complex *data_in, *data_out, *data_out2;
43c4762a1bSJed Brown ptrdiff_t alloc_local, local_n0, local_0_start;
44c4762a1bSJed Brown
45dd400576SPatrick Sanan if (rank == 0) printf("Use FFTW without PETSc-FFTW interface\n");
46c4762a1bSJed Brown fftw_mpi_init();
47c4762a1bSJed Brown N = N0 * N1;
48c4762a1bSJed Brown alloc_local = fftw_mpi_local_size_2d(N0, N1, PETSC_COMM_WORLD, &local_n0, &local_0_start);
49c4762a1bSJed Brown
50c4762a1bSJed Brown data_in = (fftw_complex *)fftw_malloc(sizeof(fftw_complex) * alloc_local);
51c4762a1bSJed Brown data_out = (fftw_complex *)fftw_malloc(sizeof(fftw_complex) * alloc_local);
52c4762a1bSJed Brown data_out2 = (fftw_complex *)fftw_malloc(sizeof(fftw_complex) * alloc_local);
53c4762a1bSJed Brown
549566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, (PetscInt)local_n0 * N1, (PetscInt)N, (const PetscScalar *)data_in, &x));
559566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x, "Real Space vector"));
569566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, (PetscInt)local_n0 * N1, (PetscInt)N, (const PetscScalar *)data_out, &y));
579566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)y, "Frequency space vector"));
589566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, (PetscInt)local_n0 * N1, (PetscInt)N, (const PetscScalar *)data_out2, &z));
599566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)z, "Reconstructed vector"));
60c4762a1bSJed Brown
61c4762a1bSJed Brown fplan = fftw_mpi_plan_dft_2d(N0, N1, data_in, data_out, PETSC_COMM_WORLD, FFTW_FORWARD, FFTW_ESTIMATE);
62c4762a1bSJed Brown bplan = fftw_mpi_plan_dft_2d(N0, N1, data_out, data_out2, PETSC_COMM_WORLD, FFTW_BACKWARD, FFTW_ESTIMATE);
63c4762a1bSJed Brown
649566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm));
659566063dSJacob Faibussowitsch if (view) PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
66c4762a1bSJed Brown
67c4762a1bSJed Brown fftw_execute(fplan);
689566063dSJacob Faibussowitsch if (view) PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
69c4762a1bSJed Brown
70c4762a1bSJed Brown fftw_execute(bplan);
71c4762a1bSJed Brown
72c4762a1bSJed Brown /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */
73c4762a1bSJed Brown a = 1.0 / (PetscReal)N;
749566063dSJacob Faibussowitsch PetscCall(VecScale(z, a));
759566063dSJacob Faibussowitsch if (view) PetscCall(VecView(z, PETSC_VIEWER_STDOUT_WORLD));
769566063dSJacob Faibussowitsch PetscCall(VecAXPY(z, -1.0, x));
779566063dSJacob Faibussowitsch PetscCall(VecNorm(z, NORM_1, &enorm));
7848a46eb9SPierre Jolivet if (enorm > 1.e-11) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Error norm of |x - z| %g\n", (double)enorm));
79c4762a1bSJed Brown
80c4762a1bSJed Brown /* Free spaces */
81c4762a1bSJed Brown fftw_destroy_plan(fplan);
82c4762a1bSJed Brown fftw_destroy_plan(bplan);
839371c9d4SSatish Balay fftw_free(data_in);
849371c9d4SSatish Balay PetscCall(VecDestroy(&x));
859371c9d4SSatish Balay fftw_free(data_out);
869371c9d4SSatish Balay PetscCall(VecDestroy(&y));
879371c9d4SSatish Balay fftw_free(data_out2);
889371c9d4SSatish Balay PetscCall(VecDestroy(&z));
89c4762a1bSJed Brown
90c4762a1bSJed Brown } else {
91c4762a1bSJed Brown /* Use PETSc-FFTW interface */
92c4762a1bSJed Brown /*-------------------------------------------*/
93c4762a1bSJed Brown PetscInt i, *dim, k, DIM;
94c4762a1bSJed Brown Mat A;
95c4762a1bSJed Brown Vec input, output;
96c4762a1bSJed Brown
97c4762a1bSJed Brown N = 30;
98c4762a1bSJed Brown for (i = 2; i < 3; i++) { /* (i=3,4: -- error in VecScatterPetscToFFTW(A,input,x); */
99c4762a1bSJed Brown DIM = i;
1009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(i, &dim));
101ad540459SPierre Jolivet for (k = 0; k < i; k++) dim[k] = 30;
102c4762a1bSJed Brown N *= dim[i - 1];
103c4762a1bSJed Brown
104c4762a1bSJed Brown /* Create FFTW object */
10548a46eb9SPierre Jolivet if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Use PETSc-FFTW interface...%d-DIM:%d \n", DIM, N));
1069566063dSJacob Faibussowitsch PetscCall(MatCreateFFT(PETSC_COMM_WORLD, DIM, dim, MATFFTW, &A));
107c4762a1bSJed Brown
108c4762a1bSJed Brown /* Create FFTW vectors that are compatible with parallel layout of A */
1099566063dSJacob Faibussowitsch PetscCall(MatCreateVecsFFTW(A, &x, &y, &z));
1109566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x, "Real space vector"));
1119566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)y, "Frequency space vector"));
1129566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)z, "Reconstructed vector"));
113c4762a1bSJed Brown
114c4762a1bSJed Brown /* Create and set PETSc vector */
1159566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &input));
1169566063dSJacob Faibussowitsch PetscCall(VecSetSizes(input, PETSC_DECIDE, N));
1179566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(input));
1189566063dSJacob Faibussowitsch PetscCall(VecSetRandom(input, rdm));
1199566063dSJacob Faibussowitsch PetscCall(VecDuplicate(input, &output));
1209566063dSJacob Faibussowitsch if (view) PetscCall(VecView(input, PETSC_VIEWER_STDOUT_WORLD));
121c4762a1bSJed Brown
122c4762a1bSJed Brown /* Vector input is copied to another vector x using VecScatterPetscToFFTW. This is because the user data
123c4762a1bSJed Brown can have any parallel layout. But FFTW requires special parallel layout of the data. Hence the original
124c4762a1bSJed Brown data which is in the vector "input" here, needs to be copied to a vector x, which has the correct parallel
125c4762a1bSJed Brown layout for FFTW. Also, during parallel real transform, this pads extra zeros automatically
126c4762a1bSJed Brown at the end of last dimension. This padding is required by FFTW to perform parallel real D.F.T. */
1279566063dSJacob Faibussowitsch PetscCall(VecScatterPetscToFFTW(A, input, x)); /* buggy for dim = 3, 4... */
128c4762a1bSJed Brown
129c4762a1bSJed Brown /* Apply FFTW_FORWARD and FFTW_BACKWARD */
1309566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, y));
1319566063dSJacob Faibussowitsch if (view) PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
1329566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, y, z));
133c4762a1bSJed Brown
134c4762a1bSJed Brown /* Output from Backward DFT needs to be modified to obtain user readable data the routine VecScatterFFTWToPetsc
135c4762a1bSJed Brown performs the job. In some sense this is the reverse operation of VecScatterPetscToFFTW. This routine gets rid of
136c4762a1bSJed Brown the extra spaces that were artificially padded to perform real parallel transform. */
1379566063dSJacob Faibussowitsch PetscCall(VecScatterFFTWToPetsc(A, z, output));
138c4762a1bSJed Brown
139c4762a1bSJed Brown /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */
140c4762a1bSJed Brown a = 1.0 / (PetscReal)N;
1419566063dSJacob Faibussowitsch PetscCall(VecScale(output, a));
1429566063dSJacob Faibussowitsch if (view) PetscCall(VecView(output, PETSC_VIEWER_STDOUT_WORLD));
1439566063dSJacob Faibussowitsch PetscCall(VecAXPY(output, -1.0, input));
1449566063dSJacob Faibussowitsch PetscCall(VecNorm(output, NORM_1, &enorm));
14548a46eb9SPierre Jolivet if (enorm > 1.e-09 && rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Error norm of |x - z| %e\n", enorm));
146c4762a1bSJed Brown
147c4762a1bSJed Brown /* Free spaces */
1489566063dSJacob Faibussowitsch PetscCall(PetscFree(dim));
1499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&input));
1509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&output));
1519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y));
1539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&z));
1549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
155c4762a1bSJed Brown }
156c4762a1bSJed Brown }
1579566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm));
1589566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
159b122ec5aSJacob Faibussowitsch return 0;
160c4762a1bSJed Brown }
161c4762a1bSJed Brown
162c4762a1bSJed Brown /*TEST
163c4762a1bSJed Brown
164c4762a1bSJed Brown build:
1650cf2e031SBarry Smith requires: !mpiuni fftw !complex
166c4762a1bSJed Brown
167c4762a1bSJed Brown test:
168c4762a1bSJed Brown output_file: output/ex158.out
169c4762a1bSJed Brown
170c4762a1bSJed Brown test:
171c4762a1bSJed Brown suffix: 2
172c4762a1bSJed Brown nsize: 3
173c4762a1bSJed Brown
174c4762a1bSJed Brown TEST*/
175