1d5b43468SJose E. Roman /* This program illustrates use of parallel real FFT*/
2c4762a1bSJed Brown static char help[] = "This program illustrates the use of parallel real 3D fftw (without PETSc interface)";
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown #include <fftw3.h>
5c4762a1bSJed Brown #include <fftw3-mpi.h>
6c4762a1bSJed Brown
main(int argc,char ** args)7d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
8d71ae5a4SJacob Faibussowitsch {
9c4762a1bSJed Brown ptrdiff_t N0 = 256, N1 = 256, N2 = 256, N3 = 2, dim[4];
10c4762a1bSJed Brown fftw_plan bplan, fplan;
11c4762a1bSJed Brown fftw_complex *out;
12c4762a1bSJed Brown double *in1, *in2;
13c4762a1bSJed Brown ptrdiff_t alloc_local, local_n0, local_0_start;
14c4762a1bSJed Brown ptrdiff_t local_n1, local_1_start;
15c4762a1bSJed Brown PetscInt i, j, indx, n1;
16c4762a1bSJed Brown PetscInt size, rank, n, N, *in, N_factor, NM;
17c4762a1bSJed Brown PetscScalar *data_fin, value1, one = 1.57, zero = 0.0;
18c4762a1bSJed Brown PetscScalar a, *x_arr, *y_arr, *z_arr, enorm;
19c4762a1bSJed Brown Vec fin, fout, fout1, ini, final;
20c4762a1bSJed Brown PetscRandom rnd;
21c4762a1bSJed Brown VecScatter vecscat, vecscat1;
22c4762a1bSJed Brown IS indx1, indx2;
23c4762a1bSJed Brown PetscInt *indx3, k, l, *indx4;
24c4762a1bSJed Brown PetscInt low, tempindx, tempindx1;
25c4762a1bSJed Brown
26327415f7SBarry Smith PetscFunctionBeginUser;
27*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
28c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
29c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "This example requires real numbers. Your current scalar type is complex");
30c4762a1bSJed Brown #endif
319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
33c4762a1bSJed Brown
34c4762a1bSJed Brown PetscRandomCreate(PETSC_COMM_WORLD, &rnd);
35c4762a1bSJed Brown
36c4762a1bSJed Brown alloc_local = fftw_mpi_local_size_3d_transposed(N0, N1, N2 / 2 + 1, 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
509371c9d4SSatish Balay N = 2 * N0 * N1 * (N2 / 2 + 1);
519371c9d4SSatish Balay N_factor = N0 * N1 * N2;
529371c9d4SSatish Balay n = 2 * local_n0 * N1 * (N2 / 2 + 1);
539371c9d4SSatish Balay n1 = local_n1 * N0 * 2 * N1;
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 VecSet(fout, zero);
68c4762a1bSJed Brown VecSet(fout1, zero);
69c4762a1bSJed Brown
70c4762a1bSJed Brown VecAssemblyBegin(fin);
71c4762a1bSJed Brown VecAssemblyEnd(fin);
72c4762a1bSJed Brown /* VecView(fin,PETSC_VIEWER_STDOUT_WORLD); */
73c4762a1bSJed Brown
74c4762a1bSJed Brown VecGetArray(fin, &x_arr);
75c4762a1bSJed Brown VecGetArray(fout1, &z_arr);
76c4762a1bSJed Brown VecGetArray(fout, &y_arr);
77c4762a1bSJed Brown
78c4762a1bSJed Brown fplan = fftw_mpi_plan_dft_r2c_3d(N0, N1, N2, (double *)x_arr, (fftw_complex *)y_arr, PETSC_COMM_WORLD, FFTW_ESTIMATE);
79c4762a1bSJed Brown bplan = fftw_mpi_plan_dft_c2r_3d(N0, N1, N2, (fftw_complex *)y_arr, (double *)z_arr, PETSC_COMM_WORLD, FFTW_ESTIMATE);
80c4762a1bSJed Brown
81c4762a1bSJed Brown fftw_execute(fplan);
82c4762a1bSJed Brown fftw_execute(bplan);
83c4762a1bSJed Brown
84c4762a1bSJed Brown VecRestoreArray(fin, &x_arr);
85c4762a1bSJed Brown VecRestoreArray(fout1, &z_arr);
86c4762a1bSJed Brown VecRestoreArray(fout, &y_arr);
87c4762a1bSJed Brown
88c4762a1bSJed Brown /* a = 1.0/(PetscReal)N_factor; */
899566063dSJacob Faibussowitsch /* PetscCall(VecScale(fout1,a)); */
90c4762a1bSJed Brown VecCreate(PETSC_COMM_WORLD, &ini);
91c4762a1bSJed Brown VecCreate(PETSC_COMM_WORLD, &final);
92c4762a1bSJed Brown VecSetSizes(ini, local_n0 * N1 * N2, N_factor);
93c4762a1bSJed Brown VecSetSizes(final, local_n0 * N1 * N2, N_factor);
94c4762a1bSJed Brown /* VecSetSizes(ini,PETSC_DECIDE,N_factor); */
95c4762a1bSJed Brown /* VecSetSizes(final,PETSC_DECIDE,N_factor); */
96c4762a1bSJed Brown VecSetFromOptions(ini);
97c4762a1bSJed Brown VecSetFromOptions(final);
98c4762a1bSJed Brown
99c4762a1bSJed Brown if (N2 % 2 == 0) NM = N2 + 2;
100c4762a1bSJed Brown else NM = N2 + 1;
101c4762a1bSJed Brown
1029566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(fin, &low, NULL));
103c4762a1bSJed Brown printf("The local index is %d from %d\n", low, rank);
1049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(local_n0 * N1 * N2, &indx3));
1059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(local_n0 * N1 * N2, &indx4));
106c4762a1bSJed Brown for (i = 0; i < local_n0; i++) {
107c4762a1bSJed Brown for (j = 0; j < N1; j++) {
108c4762a1bSJed Brown for (k = 0; k < N2; k++) {
109c4762a1bSJed Brown tempindx = i * N1 * N2 + j * N2 + k;
110c4762a1bSJed Brown tempindx1 = i * N1 * NM + j * NM + k;
111c4762a1bSJed Brown
112c4762a1bSJed Brown indx3[tempindx] = local_0_start * N1 * N2 + tempindx;
113c4762a1bSJed Brown indx4[tempindx] = low + tempindx1;
114c4762a1bSJed Brown }
115c4762a1bSJed Brown /* printf("index3 %d from proc %d is \n",indx3[tempindx],rank); */
116c4762a1bSJed Brown /* printf("index4 %d from proc %d is \n",indx4[tempindx],rank); */
117c4762a1bSJed Brown }
118c4762a1bSJed Brown }
119c4762a1bSJed Brown VecGetValues(fin, local_n0 * N1 * N2, indx4, x_arr);
120c4762a1bSJed Brown VecSetValues(ini, local_n0 * N1 * N2, indx3, x_arr, INSERT_VALUES);
121c4762a1bSJed Brown VecAssemblyBegin(ini);
122c4762a1bSJed Brown VecAssemblyEnd(ini);
123c4762a1bSJed Brown
124c4762a1bSJed Brown VecGetValues(fout1, local_n0 * N1 * N2, indx4, y_arr);
125c4762a1bSJed Brown VecSetValues(final, local_n0 * N1 * N2, indx3, y_arr, INSERT_VALUES);
126c4762a1bSJed Brown VecAssemblyBegin(final);
127c4762a1bSJed Brown VecAssemblyEnd(final);
128c4762a1bSJed Brown
129c4762a1bSJed Brown printf("The local index value is %ld from %d", local_n0 * N1 * N2, rank);
130c4762a1bSJed Brown /*
131c4762a1bSJed Brown for (i=0;i<N0;i++) {
132c4762a1bSJed Brown for (j=0;j<N1;j++) {
133c4762a1bSJed Brown indx=i*N1*NM+j*NM;
134c4762a1bSJed Brown ISCreateStride(PETSC_COMM_WORLD,N2,indx,1,&indx1);
135c4762a1bSJed Brown indx=i*N1*N2+j*N2;
136c4762a1bSJed Brown ISCreateStride(PETSC_COMM_WORLD,N2,indx,1,&indx2);
137c4762a1bSJed Brown VecScatterCreate(fin,indx1,ini,indx2,&vecscat);
138c4762a1bSJed Brown VecScatterBegin(vecscat,fin,ini,INSERT_VALUES,SCATTER_FORWARD);
139c4762a1bSJed Brown VecScatterEnd(vecscat,fin,ini,INSERT_VALUES,SCATTER_FORWARD);
140c4762a1bSJed Brown VecScatterCreate(fout1,indx1,final,indx2,&vecscat1);
141c4762a1bSJed Brown VecScatterBegin(vecscat1,fout1,final,INSERT_VALUES,SCATTER_FORWARD);
142c4762a1bSJed Brown VecScatterEnd(vecscat1,fout1,final,INSERT_VALUES,SCATTER_FORWARD);
143c4762a1bSJed Brown }
144c4762a1bSJed Brown }
145c4762a1bSJed Brown */
146c4762a1bSJed Brown a = 1.0 / (PetscReal)N_factor;
1479566063dSJacob Faibussowitsch PetscCall(VecScale(fout1, a));
1489566063dSJacob Faibussowitsch PetscCall(VecScale(final, a));
149c4762a1bSJed Brown
150c4762a1bSJed Brown VecAssemblyBegin(ini);
151c4762a1bSJed Brown VecAssemblyEnd(ini);
152c4762a1bSJed Brown
153c4762a1bSJed Brown VecAssemblyBegin(final);
154c4762a1bSJed Brown VecAssemblyEnd(final);
155c4762a1bSJed Brown
156c4762a1bSJed Brown /* VecView(final,PETSC_VIEWER_STDOUT_WORLD); */
1579566063dSJacob Faibussowitsch PetscCall(VecAXPY(final, -1.0, ini));
1589566063dSJacob Faibussowitsch PetscCall(VecNorm(final, NORM_1, &enorm));
1599566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Error norm of |x - z| = %e\n", enorm));
160c4762a1bSJed Brown fftw_destroy_plan(fplan);
161c4762a1bSJed Brown fftw_destroy_plan(bplan);
1629371c9d4SSatish Balay fftw_free(in1);
1639371c9d4SSatish Balay PetscCall(VecDestroy(&fin));
1649371c9d4SSatish Balay fftw_free(out);
1659371c9d4SSatish Balay PetscCall(VecDestroy(&fout));
1669371c9d4SSatish Balay fftw_free(in2);
1679371c9d4SSatish Balay PetscCall(VecDestroy(&fout1));
168c4762a1bSJed Brown
1699566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
170b122ec5aSJacob Faibussowitsch return 0;
171c4762a1bSJed Brown }
172