xref: /petsc/src/mat/tests/ex144.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown /* This program illustrates use of parallel real FFT */
2c4762a1bSJed Brown static char help[]="This program illustrates the use of parallel real 2D fft using fftw without PETSc interface";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown #include <fftw3.h>
6c4762a1bSJed Brown #include <fftw3-mpi.h>
7c4762a1bSJed Brown 
8c4762a1bSJed Brown int main(int argc,char **args)
9c4762a1bSJed Brown {
10c4762a1bSJed Brown   const ptrdiff_t N0=2056,N1=2056;
11c4762a1bSJed Brown   fftw_plan       bplan,fplan;
12c4762a1bSJed Brown   fftw_complex    *out;
13c4762a1bSJed Brown   double          *in1,*in2;
14c4762a1bSJed Brown   ptrdiff_t       alloc_local,local_n0,local_0_start;
15c4762a1bSJed Brown   ptrdiff_t       local_n1,local_1_start;
16c4762a1bSJed Brown   PetscInt        i,j;
17c4762a1bSJed Brown   PetscMPIInt     size,rank;
18c4762a1bSJed Brown   int             n,N,N_factor,NM;
19c4762a1bSJed Brown   PetscScalar     one=2.0,zero=0.5;
20c4762a1bSJed Brown   PetscScalar     two=4.0,three=8.0,four=16.0;
21c4762a1bSJed Brown   PetscScalar     a,*x_arr,*y_arr,*z_arr;
22c4762a1bSJed Brown   PetscReal       enorm;
23c4762a1bSJed Brown   Vec             fin,fout,fout1;
24c4762a1bSJed Brown   Vec             ini,final;
25c4762a1bSJed Brown   PetscRandom     rnd;
26c4762a1bSJed Brown   PetscInt        *indx3,tempindx,low,*indx4,tempindx1;
27c4762a1bSJed Brown 
28*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
295f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
305f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
31c4762a1bSJed Brown 
325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rnd));
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   alloc_local = fftw_mpi_local_size_2d_transposed(N0,N1/2+1,PETSC_COMM_WORLD,&local_n0,&local_0_start,&local_n1,&local_1_start);
35c4762a1bSJed Brown #if defined(DEBUGGING)
36c4762a1bSJed Brown   printf("The value alloc_local is %ld from process %d\n",alloc_local,rank);
37c4762a1bSJed Brown   printf("The value local_n0 is %ld from process %d\n",local_n0,rank);
38c4762a1bSJed Brown   printf("The value local_0_start is  %ld from process %d\n",local_0_start,rank);
39c4762a1bSJed Brown /*    printf("The value local_n1 is  %ld from process %d\n",local_n1,rank); */
40c4762a1bSJed Brown /*    printf("The value local_1_start is  %ld from process %d\n",local_1_start,rank); */
41c4762a1bSJed Brown /*    printf("The value local_n0 is  %ld from process %d\n",local_n0,rank); */
42c4762a1bSJed Brown #endif
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* Allocate space for input and output arrays  */
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/2+1);
50c4762a1bSJed Brown   N_factor = N0*N1;
51c4762a1bSJed Brown   n        = 2*local_n0*(N1/2+1);
52c4762a1bSJed Brown 
53c4762a1bSJed Brown /*    printf("The value N is  %d from process %d\n",N,rank);  */
54c4762a1bSJed Brown /*    printf("The value n is  %d from process %d\n",n,rank);  */
55c4762a1bSJed Brown /*    printf("The value n1 is  %d from process %d\n",n1,rank);*/
56c4762a1bSJed Brown   /* Creating data vector and accompanying array with VeccreateMPIWithArray */
575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,N,(PetscScalar*)in1,&fin));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,N,(PetscScalar*)out,&fout));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,N,(PetscScalar*)in2,&fout1));
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* Set the vector with random data */
625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(fin,zero));
63c4762a1bSJed Brown /*    for (i=0;i<N0*N1;i++) */
64c4762a1bSJed Brown /*       { */
65c4762a1bSJed Brown /*       VecSetValues(fin,1,&i,&one,INSERT_VALUES); */
66c4762a1bSJed Brown /*     } */
67c4762a1bSJed Brown 
68c4762a1bSJed Brown /*    VecSet(fin,one); */
69c4762a1bSJed Brown   i    =0;
705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetValues(fin,1,&i,&one,INSERT_VALUES));
71c4762a1bSJed Brown   i    =1;
725f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetValues(fin,1,&i,&two,INSERT_VALUES));
73c4762a1bSJed Brown   i    =4;
745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetValues(fin,1,&i,&three,INSERT_VALUES));
75c4762a1bSJed Brown   i    =5;
765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetValues(fin,1,&i,&four,INSERT_VALUES));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(fin));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(fin));
79c4762a1bSJed Brown 
805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(fout,zero));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(fout1,zero));
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   /* Get the meaningful portion of array */
845f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(fin,&x_arr));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(fout1,&z_arr));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(fout,&y_arr));
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   fplan=fftw_mpi_plan_dft_r2c_2d(N0,N1,(double*)x_arr,(fftw_complex*)y_arr,PETSC_COMM_WORLD,FFTW_ESTIMATE);
89c4762a1bSJed Brown   bplan=fftw_mpi_plan_dft_c2r_2d(N0,N1,(fftw_complex*)y_arr,(double*)z_arr,PETSC_COMM_WORLD,FFTW_ESTIMATE);
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   fftw_execute(fplan);
92c4762a1bSJed Brown   fftw_execute(bplan);
93c4762a1bSJed Brown 
945f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(fin,&x_arr));
955f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(fout1,&z_arr));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(fout,&y_arr));
97c4762a1bSJed Brown 
98c4762a1bSJed Brown /*    VecView(fin,PETSC_VIEWER_STDOUT_WORLD); */
995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&ini));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&final));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(ini,local_n0*N1,N0*N1));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(final,local_n0*N1,N0*N1));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(ini));
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(final));
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   if (N1%2==0) {
107c4762a1bSJed Brown     NM = N1+2;
108c4762a1bSJed Brown   } else {
109c4762a1bSJed Brown     NM = N1+1;
110c4762a1bSJed Brown   }
111c4762a1bSJed Brown   /*printf("The Value of NM is %d",NM); */
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(fin,&low,NULL));
113c4762a1bSJed Brown   /*printf("The local index is %d from %d\n",low,rank); */
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(local_n0*N1,&indx3));
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(local_n0*N1,&indx4));
116c4762a1bSJed Brown   for (i=0;i<local_n0;i++) {
117c4762a1bSJed Brown     for (j=0;j<N1;j++) {
118c4762a1bSJed Brown       tempindx  = i*N1 + j;
119c4762a1bSJed Brown       tempindx1 = i*NM + j;
120c4762a1bSJed Brown 
121c4762a1bSJed Brown       indx3[tempindx]=local_0_start*N1+tempindx;
122c4762a1bSJed Brown       indx4[tempindx]=low+tempindx1;
123c4762a1bSJed Brown       /*          printf("index3 %d from proc %d is \n",indx3[tempindx],rank); */
124c4762a1bSJed Brown       /*          printf("index4 %d from proc %d is \n",indx4[tempindx],rank); */
125c4762a1bSJed Brown     }
126c4762a1bSJed Brown   }
127c4762a1bSJed Brown 
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(local_n0*N1,&x_arr,local_n0*N1,&y_arr)); /* arr must be allocated for VecGetValues() */
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetValues(fin,local_n0*N1,indx4,(PetscScalar*)x_arr));
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetValues(ini,local_n0*N1,indx3,x_arr,INSERT_VALUES));
131c4762a1bSJed Brown 
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(ini));
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(ini));
134c4762a1bSJed Brown 
1355f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetValues(fout1,local_n0*N1,indx4,y_arr));
1365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetValues(final,local_n0*N1,indx3,y_arr,INSERT_VALUES));
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(final));
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(final));
1395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(x_arr,y_arr));
140c4762a1bSJed Brown 
141c4762a1bSJed Brown /*
142c4762a1bSJed Brown     VecScatter      vecscat;
143c4762a1bSJed Brown     IS              indx1,indx2;
144c4762a1bSJed Brown     for (i=0;i<N0;i++) {
145c4762a1bSJed Brown        indx = i*NM;
146c4762a1bSJed Brown        ISCreateStride(PETSC_COMM_WORLD,N1,indx,1,&indx1);
147c4762a1bSJed Brown        indx = i*N1;
148c4762a1bSJed Brown        ISCreateStride(PETSC_COMM_WORLD,N1,indx,1,&indx2);
149c4762a1bSJed Brown        VecScatterCreate(fin,indx1,ini,indx2,&vecscat);
150c4762a1bSJed Brown        VecScatterBegin(vecscat,fin,ini,INSERT_VALUES,SCATTER_FORWARD);
151c4762a1bSJed Brown        VecScatterEnd(vecscat,fin,ini,INSERT_VALUES,SCATTER_FORWARD);
152c4762a1bSJed Brown        VecScatterBegin(vecscat,fout1,final,INSERT_VALUES,SCATTER_FORWARD);
153c4762a1bSJed Brown        VecScatterEnd(vecscat,fout1,final,INSERT_VALUES,SCATTER_FORWARD);
154c4762a1bSJed Brown     }
155c4762a1bSJed Brown */
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   a    = 1.0/(PetscReal)N_factor;
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(fout1,a));
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(final,a));
160c4762a1bSJed Brown 
161c4762a1bSJed Brown /*    VecView(ini,PETSC_VIEWER_STDOUT_WORLD);   */
162c4762a1bSJed Brown /*    VecView(final,PETSC_VIEWER_STDOUT_WORLD); */
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(final,-1.0,ini));
164c4762a1bSJed Brown 
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(final,NORM_1,&enorm));
166c4762a1bSJed Brown   if (enorm > 1.e-10) {
1675f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"  Error norm of |x - z|  = %e\n",enorm));
168c4762a1bSJed Brown   }
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   /* Execute fftw with function fftw_execute and destroy it after execution */
171c4762a1bSJed Brown   fftw_destroy_plan(fplan);
172c4762a1bSJed Brown   fftw_destroy_plan(bplan);
1735f80ce2aSJacob Faibussowitsch   fftw_free(in1);  CHKERRQ(VecDestroy(&fin));
1745f80ce2aSJacob Faibussowitsch   fftw_free(out);  CHKERRQ(VecDestroy(&fout));
1755f80ce2aSJacob Faibussowitsch   fftw_free(in2);  CHKERRQ(VecDestroy(&fout1));
176c4762a1bSJed Brown 
1775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&ini));
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&final));
179c4762a1bSJed Brown 
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rnd));
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(indx3));
1825f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(indx4));
183*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
184*b122ec5aSJacob Faibussowitsch   return 0;
185c4762a1bSJed Brown }
186c4762a1bSJed Brown 
187c4762a1bSJed Brown /*TEST
188c4762a1bSJed Brown 
189c4762a1bSJed Brown    build:
1900cf2e031SBarry Smith       requires: !mpiuni fftw !complex
191c4762a1bSJed Brown 
192c4762a1bSJed Brown    test:
193c4762a1bSJed Brown       output_file: output/ex144.out
194c4762a1bSJed Brown 
195c4762a1bSJed Brown    test:
196c4762a1bSJed Brown       suffix: 2
197c4762a1bSJed Brown       nsize: 3
198c4762a1bSJed Brown       output_file: output/ex144.out
199c4762a1bSJed Brown 
200c4762a1bSJed Brown TEST*/
201