xref: /petsc/src/vec/is/sf/tests/ex14.c (revision dfd57a172ac9fa6c7b5fe6de6ab5df85cefc2996)
1 
2 static char help[]= "Test event log of VecScatter with various block sizes\n\n";
3 
4 #include <petscvec.h>
5 
6 int main(int argc,char **argv)
7 {
8   PetscErrorCode     ierr;
9   PetscInt           i,j,low,high,n=256,N,errors,tot_errors;
10   PetscInt           bs=1,ix[2],iy[2];
11   PetscMPIInt        nproc,rank;
12   PetscScalar       *xval;
13   const PetscScalar *yval;
14   Vec                x,y;
15   IS                 isx,isy;
16   VecScatter         ctx;
17   const PetscInt     niter = 10;
18 #if defined(PETSC_USE_LOG)
19   PetscLogStage      stage1,stage2;
20   PetscLogEvent      event1,event2;
21   PetscLogDouble     numMessages,messageLength;
22   PetscEventPerfInfo eventInfo;
23   PetscInt           tot_msg,tot_len,avg_len;
24 #endif
25 
26   PetscFunctionBegin;
27   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
28   ierr = PetscLogDefaultBegin();CHKERRQ(ierr);
29   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&nproc);CHKERRMPI(ierr);
30   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
31 
32   ierr = PetscLogStageRegister("Scatter(bs=1)", &stage1);CHKERRQ(ierr);
33   ierr = PetscLogEventRegister("VecScatter(bs=1)", PETSC_OBJECT_CLASSID, &event1);CHKERRQ(ierr);
34   ierr = PetscLogStageRegister("Scatter(bs=4)", &stage2);CHKERRQ(ierr);
35   ierr = PetscLogEventRegister("VecScatter(bs=4)", PETSC_OBJECT_CLASSID, &event2);CHKERRQ(ierr);
36 
37   /* Create a parallel vector x and a sequential vector y */
38   ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
39   ierr = VecSetSizes(x,n,PETSC_DECIDE);CHKERRQ(ierr);
40   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
41   ierr = VecGetOwnershipRange(x,&low,&high);CHKERRQ(ierr);
42   ierr = VecGetSize(x,&N);CHKERRQ(ierr);
43   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&y);CHKERRQ(ierr);
44 
45   /*=======================================
46      test VecScatter with bs = 1
47     ======================================*/
48 
49   /* the code works as if we are going to do 3-point stencil computations on a 1D domain x,
50      which has periodic boundary conditions but the two ghosts are scatterred to beginning of y.
51    */
52   bs    = 1;
53   ix[0] = rank? low-1 : N-1; /* ix[] contains global indices of the two ghost points */
54   ix[1] = (rank != nproc-1)? high : 0;
55   iy[0] = 0;
56   iy[1] = 1;
57 
58   ierr = ISCreateGeneral(PETSC_COMM_SELF,2,ix,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr);
59   ierr = ISCreateGeneral(PETSC_COMM_SELF,2,iy,PETSC_COPY_VALUES,&isy);CHKERRQ(ierr);
60   ierr = VecScatterCreate(x,isx,y,isy,&ctx);CHKERRQ(ierr);
61   ierr = VecScatterSetUp(ctx);CHKERRQ(ierr);
62 
63   ierr = PetscLogStagePush(stage1);CHKERRQ(ierr);
64   ierr = PetscLogEventBegin(event1,0,0,0,0);CHKERRQ(ierr);
65   errors = 0;
66   for (i=0; i<niter; i++) {
67     /* set x = 0+i, 1+i, 2+i, ..., N-1+i */
68     ierr = VecGetArray(x,&xval);CHKERRQ(ierr);
69     for (j=0; j<n; j++) xval[j] = (PetscScalar)(low+j+i);
70     ierr = VecRestoreArray(x,&xval);CHKERRQ(ierr);
71     /* scatter the ghosts to y */
72     ierr = VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
73     ierr = VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
74     /* check if y has correct values */
75     ierr = VecGetArrayRead(y,&yval);CHKERRQ(ierr);
76     if ((PetscInt)PetscRealPart(yval[0]) != ix[0]+i) errors++;
77     if ((PetscInt)PetscRealPart(yval[1]) != ix[1]+i) errors++;
78     ierr = VecRestoreArrayRead(y,&yval);CHKERRQ(ierr);
79   }
80   ierr = PetscLogEventEnd(event1,0,0,0,0);CHKERRQ(ierr);
81   ierr = PetscLogStagePop();CHKERRQ(ierr);
82 
83   /* check if we found wrong values on any processors */
84   ierr = MPI_Allreduce(&errors,&tot_errors,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr);
85   if (tot_errors) { ierr = PetscPrintf(PETSC_COMM_WORLD,"Error: wrong values were scatterred in vecscatter with bs = %d\n",bs);CHKERRQ(ierr); }
86 
87   /* print out event log of VecScatter(bs=1) */
88 #if defined(PETSC_USE_LOG)
89   ierr    = PetscLogEventGetPerfInfo(stage1,event1,&eventInfo);CHKERRQ(ierr);
90   ierr    = MPI_Allreduce(&eventInfo.numMessages,  &numMessages,  1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr);
91   ierr    = MPI_Allreduce(&eventInfo.messageLength,&messageLength,1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr);
92   tot_msg = (PetscInt)numMessages*0.5; /* two MPI calls (Send & Recv) per message */
93   tot_len = (PetscInt)messageLength*0.5;
94   avg_len = tot_msg? (PetscInt)(messageLength/numMessages) : 0;
95   /* when nproc > 2, tot_msg = 2*nproc*niter, tot_len = tot_msg*sizeof(PetscScalar)*bs */
96   ierr    = PetscPrintf(PETSC_COMM_WORLD,"VecScatter(bs=%d) has sent out %d messages, total %d bytes, with average length %d bytes\n",bs,tot_msg,tot_len,avg_len);CHKERRQ(ierr);
97 #endif
98 
99   ierr = ISDestroy(&isx);CHKERRQ(ierr);
100   ierr = ISDestroy(&isy);CHKERRQ(ierr);
101   ierr = VecScatterDestroy(&ctx);CHKERRQ(ierr);
102 
103   /*=======================================
104      test VecScatter with bs = 4
105     ======================================*/
106 
107   /* similar to the 3-point stencil above, except that this time a ghost is a block */
108   bs    = 4; /* n needs to be a multiple of bs to make the following code work */
109   ix[0] = rank? low/bs-1 : N/bs-1; /* ix[] contains global indices of the two ghost blocks */
110   ix[1] = (rank != nproc-1)? high/bs : 0;
111   iy[0] = 0;
112   iy[1] = 1;
113 
114   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2,ix,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr);
115   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2,iy,PETSC_COPY_VALUES,&isy);CHKERRQ(ierr);
116 
117   ierr = VecScatterCreate(x,isx,y,isy,&ctx);CHKERRQ(ierr);
118    /* Call SetUp explicitly, otherwise messages in implicit SetUp will be counted in events below */
119   ierr = VecScatterSetUp(ctx);CHKERRQ(ierr);
120 
121   ierr = PetscLogStagePush(stage2);CHKERRQ(ierr);
122   ierr = PetscLogEventBegin(event2,0,0,0,0);CHKERRQ(ierr);
123   errors = 0;
124   for (i=0; i<niter; i++) {
125     /* set x = 0+i, 1+i, 2+i, ..., N-1+i */
126     ierr = VecGetArray(x,&xval);CHKERRQ(ierr);
127     for (j=0; j<n; j++) xval[j] = (PetscScalar)(low+j+i);
128     ierr = VecRestoreArray(x,&xval);CHKERRQ(ierr);
129     /* scatter the ghost blocks to y */
130     ierr = VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
131     ierr = VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
132     /* check if y has correct values */
133     ierr = VecGetArrayRead(y,&yval);CHKERRQ(ierr);
134     if ((PetscInt)PetscRealPart(yval[0])  != ix[0]*bs+i) errors++;
135     if ((PetscInt)PetscRealPart(yval[bs]) != ix[1]*bs+i) errors++;
136     ierr = VecRestoreArrayRead(y,&yval);CHKERRQ(ierr);
137   }
138   ierr = PetscLogEventEnd(event2,0,0,0,0);CHKERRQ(ierr);
139   ierr = PetscLogStagePop();CHKERRQ(ierr);
140 
141   /* check if we found wrong values on any processors */
142   ierr = MPI_Allreduce(&errors,&tot_errors,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr);
143   if (tot_errors) { ierr = PetscPrintf(PETSC_COMM_WORLD,"Error: wrong values were scatterred in vecscatter with bs = %d\n",bs);CHKERRQ(ierr); }
144 
145   /* print out event log of VecScatter(bs=4) */
146 #if defined(PETSC_USE_LOG)
147   ierr    = PetscLogEventGetPerfInfo(stage2,event2,&eventInfo);CHKERRQ(ierr);
148   ierr    = MPI_Allreduce(&eventInfo.numMessages,  &numMessages,  1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr);
149   ierr    = MPI_Allreduce(&eventInfo.messageLength,&messageLength,1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr);
150   tot_msg = (PetscInt)numMessages*0.5; /* two MPI calls (Send & Recv) per message */
151   tot_len = (PetscInt)messageLength*0.5;
152   avg_len = tot_msg? (PetscInt)(messageLength/numMessages) : 0;
153   /* when nproc > 2, tot_msg = 2*nproc*niter, tot_len = tot_msg*sizeof(PetscScalar)*bs */
154   ierr    = PetscPrintf(PETSC_COMM_WORLD,"VecScatter(bs=%d) has sent out %d messages, total %d bytes, with average length %d bytes\n",bs,tot_msg,tot_len,avg_len);CHKERRQ(ierr);
155 #endif
156 
157   ierr = PetscPrintf(PETSC_COMM_WORLD,"Program finished\n");CHKERRQ(ierr);
158   ierr = ISDestroy(&isx);CHKERRQ(ierr);
159   ierr = ISDestroy(&isy);CHKERRQ(ierr);
160   ierr = VecScatterDestroy(&ctx);CHKERRQ(ierr);
161   ierr = VecDestroy(&x);CHKERRQ(ierr);
162   ierr = VecDestroy(&y);CHKERRQ(ierr);
163   ierr = PetscFinalize();
164   return ierr;
165 }
166 
167 /*TEST
168 
169    test:
170       nsize: 4
171       args:
172       requires: double defined(PETSC_USE_LOG)
173 
174 TEST*/
175