xref: /petsc/src/mat/tests/ex147.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown /* This program illustrates use of parallel real FFT */
2c4762a1bSJed Brown static char help[] = "This program illustrates the use of parallel real multi-dimensional fftw (without PETSc interface)";
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown #include <fftw3.h>
5c4762a1bSJed Brown #include <fftw3-mpi.h>
6c4762a1bSJed Brown 
7*9371c9d4SSatish Balay int main(int argc, char **args) {
8c4762a1bSJed Brown   ptrdiff_t     N0 = 2, N1 = 2, N2 = 2, N3 = 2, dim[4], N, D;
9c4762a1bSJed Brown   fftw_plan     bplan, fplan;
10c4762a1bSJed Brown   fftw_complex *out;
11c4762a1bSJed Brown   double       *in1, *in2;
12c4762a1bSJed Brown   ptrdiff_t     alloc_local, local_n0, local_0_start;
13c4762a1bSJed Brown   ptrdiff_t     local_n1, local_1_start;
14c4762a1bSJed Brown   PetscInt      i, j, indx[100], n1;
15c4762a1bSJed Brown   PetscInt      size, rank, n, *in, N_factor;
16c4762a1bSJed Brown   PetscScalar  *data_fin, value1, one = 1.0, zero = 0.0;
17c4762a1bSJed Brown   PetscScalar   a, *x_arr, *y_arr, *z_arr, enorm;
18c4762a1bSJed Brown   Vec           fin, fout, fout1, x, y;
19c4762a1bSJed Brown   PetscRandom   rnd;
20c4762a1bSJed Brown 
21327415f7SBarry Smith   PetscFunctionBeginUser;
229566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
23c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
24c4762a1bSJed Brown   SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "This example requires real numbers. Your current scalar type is complex");
25c4762a1bSJed Brown #endif
269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   PetscRandomCreate(PETSC_COMM_WORLD, &rnd);
30c4762a1bSJed Brown   D      = 4;
31*9371c9d4SSatish Balay   dim[0] = N0;
32*9371c9d4SSatish Balay   dim[1] = N1;
33*9371c9d4SSatish Balay   dim[2] = N2;
34*9371c9d4SSatish Balay   dim[3] = N3 / 2 + 1;
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   alloc_local = fftw_mpi_local_size_transposed(D, dim, PETSC_COMM_WORLD, &local_n0, &local_0_start, &local_n1, &local_1_start);
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   printf("The value alloc_local is %ld from process %d\n", alloc_local, rank);
39c4762a1bSJed Brown   printf("The value local_n0 is %ld from process %d\n", local_n0, rank);
40c4762a1bSJed Brown   printf("The value local_0_start is  %ld from process %d\n", local_0_start, rank);
41c4762a1bSJed Brown   printf("The value local_n1 is  %ld from process %d\n", local_n1, rank);
42c4762a1bSJed Brown   printf("The value local_1_start is  %ld from process %d\n", local_1_start, rank);
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* Allocate space for input and output arrays  */
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   in1 = (double *)fftw_malloc(sizeof(double) * alloc_local * 2);
47c4762a1bSJed Brown   in2 = (double *)fftw_malloc(sizeof(double) * alloc_local * 2);
48c4762a1bSJed Brown   out = (fftw_complex *)fftw_malloc(sizeof(fftw_complex) * alloc_local);
49c4762a1bSJed Brown 
50*9371c9d4SSatish Balay   N        = 2 * N0 * N1 * N2 * (N3 / 2 + 1);
51*9371c9d4SSatish Balay   N_factor = N0 * N1 * N2 * N3;
52*9371c9d4SSatish Balay   n        = 2 * local_n0 * N1 * N2 * (N3 / 2 + 1);
53*9371c9d4SSatish Balay   n1       = local_n1 * N0 * 2 * N1 * N2;
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   /*    printf("The value N is  %d from process %d\n",N,rank); */
56c4762a1bSJed Brown   /*    printf("The value n is  %d from process %d\n",n,rank); */
57c4762a1bSJed Brown   /*    printf("The value n1 is  %d from process %d\n",n1,rank); */
58c4762a1bSJed Brown   /* Creating data vector and accompanying array with VeccreateMPIWithArray */
599566063dSJacob Faibussowitsch   PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, n, N, (PetscScalar *)in1, &fin));
609566063dSJacob Faibussowitsch   PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, n, N, (PetscScalar *)out, &fout));
619566063dSJacob Faibussowitsch   PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, n, N, (PetscScalar *)in2, &fout1));
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   /*    VecGetSize(fin,&size); */
64c4762a1bSJed Brown   /*    printf("The size is %d\n",size); */
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   VecSet(fin, one);
67c4762a1bSJed Brown   /*    VecAssemblyBegin(fin); */
68c4762a1bSJed Brown   /*    VecAssemblyEnd(fin); */
69c4762a1bSJed Brown   /*    VecView(fin,PETSC_VIEWER_STDOUT_WORLD); */
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   VecGetArray(fin, &x_arr);
72c4762a1bSJed Brown   VecGetArray(fout1, &z_arr);
73c4762a1bSJed Brown   VecGetArray(fout, &y_arr);
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   dim[3] = N3;
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   fplan = fftw_mpi_plan_dft_r2c(D, dim, (double *)x_arr, (fftw_complex *)y_arr, PETSC_COMM_WORLD, FFTW_ESTIMATE);
78c4762a1bSJed Brown   bplan = fftw_mpi_plan_dft_c2r(D, dim, (fftw_complex *)y_arr, (double *)z_arr, PETSC_COMM_WORLD, FFTW_ESTIMATE);
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   fftw_execute(fplan);
81c4762a1bSJed Brown   fftw_execute(bplan);
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   VecRestoreArray(fin, &x_arr);
84c4762a1bSJed Brown   VecRestoreArray(fout1, &z_arr);
85c4762a1bSJed Brown   VecRestoreArray(fout, &y_arr);
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /*    a = 1.0/(PetscReal)N_factor; */
889566063dSJacob Faibussowitsch   /*    PetscCall(VecScale(fout1,a)); */
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   VecAssemblyBegin(fout1);
91c4762a1bSJed Brown   VecAssemblyEnd(fout1);
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   VecView(fout1, PETSC_VIEWER_STDOUT_WORLD);
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   fftw_destroy_plan(fplan);
96c4762a1bSJed Brown   fftw_destroy_plan(bplan);
97*9371c9d4SSatish Balay   fftw_free(in1);
98*9371c9d4SSatish Balay   PetscCall(VecDestroy(&fin));
99*9371c9d4SSatish Balay   fftw_free(out);
100*9371c9d4SSatish Balay   PetscCall(VecDestroy(&fout));
101*9371c9d4SSatish Balay   fftw_free(in2);
102*9371c9d4SSatish Balay   PetscCall(VecDestroy(&fout1));
103c4762a1bSJed Brown 
1049566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
105b122ec5aSJacob Faibussowitsch   return 0;
106c4762a1bSJed Brown }
107