xref: /petsc/src/mat/tests/ex146.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown /* This program illustrates use of paralllel 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 
7ad0e0ea3SStefano Zampini int main(int argc,char **args)
8c4762a1bSJed Brown {
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 
26*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
27c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
28c4762a1bSJed Brown   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This example requires real numbers. Your current scalar type is complex");
29c4762a1bSJed Brown #endif
305f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
315f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   PetscRandomCreate(PETSC_COMM_WORLD,&rnd);
34c4762a1bSJed Brown 
35c4762a1bSJed 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);
36c4762a1bSJed Brown 
37c4762a1bSJed Brown /*    printf("The value alloc_local is %ld from process %d\n",alloc_local,rank);     */
38c4762a1bSJed Brown   printf("The value local_n0 is %ld from process %d\n",local_n0,rank);
39c4762a1bSJed Brown /*    printf("The value local_0_start is  %ld from process %d\n",local_0_start,rank);*/
40c4762a1bSJed Brown /*    printf("The value local_n1 is  %ld from process %d\n",local_n1,rank);          */
41c4762a1bSJed Brown /*    printf("The value local_1_start is  %ld from process %d\n",local_1_start,rank);*/
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   /* Allocate space for input and output arrays  */
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   in1=(double*)fftw_malloc(sizeof(double)*alloc_local*2);
46c4762a1bSJed Brown   in2=(double*)fftw_malloc(sizeof(double)*alloc_local*2);
47c4762a1bSJed Brown   out=(fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   N=2*N0*N1*(N2/2+1);N_factor=N0*N1*N2;
50c4762a1bSJed Brown   n=2*local_n0*N1*(N2/2+1);n1=local_n1*N0*2*N1;
51c4762a1bSJed Brown 
52c4762a1bSJed Brown /*    printf("The value N is  %d from process %d\n",N,rank);   */
53c4762a1bSJed Brown /*    printf("The value n is  %d from process %d\n",n,rank);   */
54c4762a1bSJed Brown /*    printf("The value n1 is  %d from process %d\n",n1,rank); */
55c4762a1bSJed Brown   /* Creating data vector and accompanying array with VeccreateMPIWithArray */
565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,N,(PetscScalar*)in1,&fin));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,N,(PetscScalar*)out,&fout));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,N,(PetscScalar*)in2,&fout1));
59c4762a1bSJed Brown 
60c4762a1bSJed Brown /*    VecGetSize(fin,&size); */
61c4762a1bSJed Brown /*    printf("The size is %d\n",size); */
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   VecSet(fin,one);
64c4762a1bSJed Brown   VecSet(fout,zero);
65c4762a1bSJed Brown   VecSet(fout1,zero);
66c4762a1bSJed Brown 
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   fplan=fftw_mpi_plan_dft_r2c_3d(N0,N1,N2,(double*)x_arr,(fftw_complex*)y_arr,PETSC_COMM_WORLD,FFTW_ESTIMATE);
76c4762a1bSJed Brown   bplan=fftw_mpi_plan_dft_c2r_3d(N0,N1,N2,(fftw_complex*)y_arr,(double*)z_arr,PETSC_COMM_WORLD,FFTW_ESTIMATE);
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   fftw_execute(fplan);
79c4762a1bSJed Brown   fftw_execute(bplan);
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   VecRestoreArray(fin,&x_arr);
82c4762a1bSJed Brown   VecRestoreArray(fout1,&z_arr);
83c4762a1bSJed Brown   VecRestoreArray(fout,&y_arr);
84c4762a1bSJed Brown 
85c4762a1bSJed Brown /*    a = 1.0/(PetscReal)N_factor; */
865f80ce2aSJacob Faibussowitsch /*    CHKERRQ(VecScale(fout1,a)); */
87c4762a1bSJed Brown   VecCreate(PETSC_COMM_WORLD,&ini);
88c4762a1bSJed Brown   VecCreate(PETSC_COMM_WORLD,&final);
89c4762a1bSJed Brown   VecSetSizes(ini,local_n0*N1*N2,N_factor);
90c4762a1bSJed Brown   VecSetSizes(final,local_n0*N1*N2,N_factor);
91c4762a1bSJed Brown /*    VecSetSizes(ini,PETSC_DECIDE,N_factor); */
92c4762a1bSJed Brown /*    VecSetSizes(final,PETSC_DECIDE,N_factor); */
93c4762a1bSJed Brown   VecSetFromOptions(ini);
94c4762a1bSJed Brown   VecSetFromOptions(final);
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   if (N2%2==0) NM=N2+2;
97c4762a1bSJed Brown   else NM=N2+1;
98c4762a1bSJed Brown 
995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(fin,&low,NULL));
100c4762a1bSJed Brown   printf("The local index is %d from %d\n",low,rank);
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(local_n0*N1*N2,&indx3));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(local_n0*N1*N2,&indx4));
103c4762a1bSJed Brown   for (i=0; i<local_n0; i++) {
104c4762a1bSJed Brown     for (j=0;j<N1;j++) {
105c4762a1bSJed Brown       for (k=0;k<N2;k++) {
106c4762a1bSJed Brown         tempindx  = i*N1*N2 + j*N2 + k;
107c4762a1bSJed Brown         tempindx1 = i*N1*NM + j*NM + k;
108c4762a1bSJed Brown 
109c4762a1bSJed Brown         indx3[tempindx]=local_0_start*N1*N2+tempindx;
110c4762a1bSJed Brown         indx4[tempindx]=low+tempindx1;
111c4762a1bSJed Brown       }
112c4762a1bSJed Brown       /*          printf("index3 %d from proc %d is \n",indx3[tempindx],rank); */
113c4762a1bSJed Brown       /*          printf("index4 %d from proc %d is \n",indx4[tempindx],rank); */
114c4762a1bSJed Brown     }
115c4762a1bSJed Brown   }
116c4762a1bSJed Brown   VecGetValues(fin,local_n0*N1*N2,indx4,x_arr);
117c4762a1bSJed Brown   VecSetValues(ini,local_n0*N1*N2,indx3,x_arr,INSERT_VALUES);
118c4762a1bSJed Brown   VecAssemblyBegin(ini);
119c4762a1bSJed Brown   VecAssemblyEnd(ini);
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   VecGetValues(fout1,local_n0*N1*N2,indx4,y_arr);
122c4762a1bSJed Brown   VecSetValues(final,local_n0*N1*N2,indx3,y_arr,INSERT_VALUES);
123c4762a1bSJed Brown   VecAssemblyBegin(final);
124c4762a1bSJed Brown   VecAssemblyEnd(final);
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   printf("The local index value is %ld from %d",local_n0*N1*N2,rank);
127c4762a1bSJed Brown /*
128c4762a1bSJed Brown   for (i=0;i<N0;i++) {
129c4762a1bSJed Brown      for (j=0;j<N1;j++) {
130c4762a1bSJed Brown         indx=i*N1*NM+j*NM;
131c4762a1bSJed Brown         ISCreateStride(PETSC_COMM_WORLD,N2,indx,1,&indx1);
132c4762a1bSJed Brown         indx=i*N1*N2+j*N2;
133c4762a1bSJed Brown         ISCreateStride(PETSC_COMM_WORLD,N2,indx,1,&indx2);
134c4762a1bSJed Brown         VecScatterCreate(fin,indx1,ini,indx2,&vecscat);
135c4762a1bSJed Brown         VecScatterBegin(vecscat,fin,ini,INSERT_VALUES,SCATTER_FORWARD);
136c4762a1bSJed Brown         VecScatterEnd(vecscat,fin,ini,INSERT_VALUES,SCATTER_FORWARD);
137c4762a1bSJed Brown         VecScatterCreate(fout1,indx1,final,indx2,&vecscat1);
138c4762a1bSJed Brown         VecScatterBegin(vecscat1,fout1,final,INSERT_VALUES,SCATTER_FORWARD);
139c4762a1bSJed Brown         VecScatterEnd(vecscat1,fout1,final,INSERT_VALUES,SCATTER_FORWARD);
140c4762a1bSJed Brown      }
141c4762a1bSJed Brown   }
142c4762a1bSJed Brown */
143c4762a1bSJed Brown   a    = 1.0/(PetscReal)N_factor;
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(fout1,a));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(final,a));
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   VecAssemblyBegin(ini);
148c4762a1bSJed Brown   VecAssemblyEnd(ini);
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   VecAssemblyBegin(final);
151c4762a1bSJed Brown   VecAssemblyEnd(final);
152c4762a1bSJed Brown 
153c4762a1bSJed Brown /*    VecView(final,PETSC_VIEWER_STDOUT_WORLD); */
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(final,-1.0,ini));
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(final,NORM_1,&enorm));
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"  Error norm of |x - z|  = %e\n",enorm));
157c4762a1bSJed Brown   fftw_destroy_plan(fplan);
158c4762a1bSJed Brown   fftw_destroy_plan(bplan);
1595f80ce2aSJacob Faibussowitsch   fftw_free(in1); CHKERRQ(VecDestroy(&fin));
1605f80ce2aSJacob Faibussowitsch   fftw_free(out); CHKERRQ(VecDestroy(&fout));
1615f80ce2aSJacob Faibussowitsch   fftw_free(in2); CHKERRQ(VecDestroy(&fout1));
162c4762a1bSJed Brown 
163*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
164*b122ec5aSJacob Faibussowitsch   return 0;
165c4762a1bSJed Brown }
166