xref: /petsc/src/vec/is/sf/tests/ex14.c (revision 9566063d113dddea24716c546802770db7481bc0)
197929ea7SJunchao Zhang 
297929ea7SJunchao Zhang static char help[]= "Test event log of VecScatter with various block sizes\n\n";
397929ea7SJunchao Zhang 
497929ea7SJunchao Zhang #include <petscvec.h>
597929ea7SJunchao Zhang 
697929ea7SJunchao Zhang int main(int argc,char **argv)
797929ea7SJunchao Zhang {
897929ea7SJunchao Zhang   PetscInt           i,j,low,high,n=256,N,errors,tot_errors;
997929ea7SJunchao Zhang   PetscInt           bs=1,ix[2],iy[2];
1097929ea7SJunchao Zhang   PetscMPIInt        nproc,rank;
1197929ea7SJunchao Zhang   PetscScalar       *xval;
1297929ea7SJunchao Zhang   const PetscScalar *yval;
1397929ea7SJunchao Zhang   Vec                x,y;
1497929ea7SJunchao Zhang   IS                 isx,isy;
1597929ea7SJunchao Zhang   VecScatter         ctx;
1697929ea7SJunchao Zhang   const PetscInt     niter = 10;
1797929ea7SJunchao Zhang #if defined(PETSC_USE_LOG)
1897929ea7SJunchao Zhang   PetscLogStage      stage1,stage2;
1997929ea7SJunchao Zhang   PetscLogEvent      event1,event2;
2097929ea7SJunchao Zhang   PetscLogDouble     numMessages,messageLength;
2197929ea7SJunchao Zhang   PetscEventPerfInfo eventInfo;
2297929ea7SJunchao Zhang   PetscInt           tot_msg,tot_len,avg_len;
2397929ea7SJunchao Zhang #endif
2497929ea7SJunchao Zhang 
2597929ea7SJunchao Zhang   PetscFunctionBegin;
26*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
27*9566063dSJacob Faibussowitsch   PetscCall(PetscLogDefaultBegin());
28*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&nproc));
29*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
3097929ea7SJunchao Zhang 
31*9566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("Scatter(bs=1)", &stage1));
32*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("VecScatter(bs=1)", PETSC_OBJECT_CLASSID, &event1));
33*9566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("Scatter(bs=4)", &stage2));
34*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("VecScatter(bs=4)", PETSC_OBJECT_CLASSID, &event2));
3597929ea7SJunchao Zhang 
3697929ea7SJunchao Zhang   /* Create a parallel vector x and a sequential vector y */
37*9566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
38*9566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x,n,PETSC_DECIDE));
39*9566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
40*9566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x,&low,&high));
41*9566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x,&N));
42*9566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&y));
4397929ea7SJunchao Zhang 
4497929ea7SJunchao Zhang   /*=======================================
4597929ea7SJunchao Zhang      test VecScatter with bs = 1
4697929ea7SJunchao Zhang     ======================================*/
4797929ea7SJunchao Zhang 
4897929ea7SJunchao Zhang   /* the code works as if we are going to do 3-point stencil computations on a 1D domain x,
4997929ea7SJunchao Zhang      which has periodic boundary conditions but the two ghosts are scatterred to beginning of y.
5097929ea7SJunchao Zhang    */
5197929ea7SJunchao Zhang   bs    = 1;
5297929ea7SJunchao Zhang   ix[0] = rank? low-1 : N-1; /* ix[] contains global indices of the two ghost points */
5397929ea7SJunchao Zhang   ix[1] = (rank != nproc-1)? high : 0;
5497929ea7SJunchao Zhang   iy[0] = 0;
5597929ea7SJunchao Zhang   iy[1] = 1;
5697929ea7SJunchao Zhang 
57*9566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,2,ix,PETSC_COPY_VALUES,&isx));
58*9566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,2,iy,PETSC_COPY_VALUES,&isy));
59*9566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x,isx,y,isy,&ctx));
60*9566063dSJacob Faibussowitsch   PetscCall(VecScatterSetUp(ctx));
6197929ea7SJunchao Zhang 
62*9566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage1));
63*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(event1,0,0,0,0));
6497929ea7SJunchao Zhang   errors = 0;
6597929ea7SJunchao Zhang   for (i=0; i<niter; i++) {
6697929ea7SJunchao Zhang     /* set x = 0+i, 1+i, 2+i, ..., N-1+i */
67*9566063dSJacob Faibussowitsch     PetscCall(VecGetArray(x,&xval));
6897929ea7SJunchao Zhang     for (j=0; j<n; j++) xval[j] = (PetscScalar)(low+j+i);
69*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(x,&xval));
7097929ea7SJunchao Zhang     /* scatter the ghosts to y */
71*9566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD));
72*9566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD));
7397929ea7SJunchao Zhang     /* check if y has correct values */
74*9566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(y,&yval));
7597929ea7SJunchao Zhang     if ((PetscInt)PetscRealPart(yval[0]) != ix[0]+i) errors++;
7697929ea7SJunchao Zhang     if ((PetscInt)PetscRealPart(yval[1]) != ix[1]+i) errors++;
77*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(y,&yval));
7897929ea7SJunchao Zhang   }
79*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(event1,0,0,0,0));
80*9566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
8197929ea7SJunchao Zhang 
8297929ea7SJunchao Zhang   /* check if we found wrong values on any processors */
83*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&errors,&tot_errors,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD));
84*9566063dSJacob Faibussowitsch   if (tot_errors) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error: wrong values were scatterred in vecscatter with bs = %" PetscInt_FMT "\n",bs));
8597929ea7SJunchao Zhang 
8697929ea7SJunchao Zhang   /* print out event log of VecScatter(bs=1) */
8797929ea7SJunchao Zhang #if defined(PETSC_USE_LOG)
88*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventGetPerfInfo(stage1,event1,&eventInfo));
89*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&eventInfo.numMessages,  &numMessages,  1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD));
90*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&eventInfo.messageLength,&messageLength,1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD));
9197929ea7SJunchao Zhang   tot_msg = (PetscInt)numMessages*0.5; /* two MPI calls (Send & Recv) per message */
9297929ea7SJunchao Zhang   tot_len = (PetscInt)messageLength*0.5;
9397929ea7SJunchao Zhang   avg_len = tot_msg? (PetscInt)(messageLength/numMessages) : 0;
9497929ea7SJunchao Zhang   /* when nproc > 2, tot_msg = 2*nproc*niter, tot_len = tot_msg*sizeof(PetscScalar)*bs */
95*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"VecScatter(bs=%" PetscInt_FMT ") has sent out %" PetscInt_FMT " messages, total %" PetscInt_FMT " bytes, with average length %" PetscInt_FMT " bytes\n",bs,tot_msg,tot_len,avg_len));
9697929ea7SJunchao Zhang #endif
9797929ea7SJunchao Zhang 
98*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isx));
99*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isy));
100*9566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&ctx));
10197929ea7SJunchao Zhang 
10297929ea7SJunchao Zhang   /*=======================================
10397929ea7SJunchao Zhang      test VecScatter with bs = 4
10497929ea7SJunchao Zhang     ======================================*/
10597929ea7SJunchao Zhang 
10697929ea7SJunchao Zhang   /* similar to the 3-point stencil above, except that this time a ghost is a block */
10797929ea7SJunchao Zhang   bs    = 4; /* n needs to be a multiple of bs to make the following code work */
10897929ea7SJunchao Zhang   ix[0] = rank? low/bs-1 : N/bs-1; /* ix[] contains global indices of the two ghost blocks */
10997929ea7SJunchao Zhang   ix[1] = (rank != nproc-1)? high/bs : 0;
11097929ea7SJunchao Zhang   iy[0] = 0;
11197929ea7SJunchao Zhang   iy[1] = 1;
11297929ea7SJunchao Zhang 
113*9566063dSJacob Faibussowitsch   PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,2,ix,PETSC_COPY_VALUES,&isx));
114*9566063dSJacob Faibussowitsch   PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,2,iy,PETSC_COPY_VALUES,&isy));
11597929ea7SJunchao Zhang 
116*9566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x,isx,y,isy,&ctx));
11797929ea7SJunchao Zhang    /* Call SetUp explicitly, otherwise messages in implicit SetUp will be counted in events below */
118*9566063dSJacob Faibussowitsch   PetscCall(VecScatterSetUp(ctx));
11997929ea7SJunchao Zhang 
120*9566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage2));
121*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(event2,0,0,0,0));
12297929ea7SJunchao Zhang   errors = 0;
12397929ea7SJunchao Zhang   for (i=0; i<niter; i++) {
12497929ea7SJunchao Zhang     /* set x = 0+i, 1+i, 2+i, ..., N-1+i */
125*9566063dSJacob Faibussowitsch     PetscCall(VecGetArray(x,&xval));
12697929ea7SJunchao Zhang     for (j=0; j<n; j++) xval[j] = (PetscScalar)(low+j+i);
127*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(x,&xval));
12897929ea7SJunchao Zhang     /* scatter the ghost blocks to y */
129*9566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD));
130*9566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD));
13197929ea7SJunchao Zhang     /* check if y has correct values */
132*9566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(y,&yval));
13397929ea7SJunchao Zhang     if ((PetscInt)PetscRealPart(yval[0])  != ix[0]*bs+i) errors++;
13497929ea7SJunchao Zhang     if ((PetscInt)PetscRealPart(yval[bs]) != ix[1]*bs+i) errors++;
135*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(y,&yval));
13697929ea7SJunchao Zhang   }
137*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(event2,0,0,0,0));
138*9566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
13997929ea7SJunchao Zhang 
14097929ea7SJunchao Zhang   /* check if we found wrong values on any processors */
141*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&errors,&tot_errors,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD));
142*9566063dSJacob Faibussowitsch   if (tot_errors) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error: wrong values were scatterred in vecscatter with bs = %" PetscInt_FMT "\n",bs));
14397929ea7SJunchao Zhang 
14497929ea7SJunchao Zhang   /* print out event log of VecScatter(bs=4) */
14597929ea7SJunchao Zhang #if defined(PETSC_USE_LOG)
146*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventGetPerfInfo(stage2,event2,&eventInfo));
147*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&eventInfo.numMessages,  &numMessages,  1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD));
148*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&eventInfo.messageLength,&messageLength,1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD));
14997929ea7SJunchao Zhang   tot_msg = (PetscInt)numMessages*0.5; /* two MPI calls (Send & Recv) per message */
15097929ea7SJunchao Zhang   tot_len = (PetscInt)messageLength*0.5;
15197929ea7SJunchao Zhang   avg_len = tot_msg? (PetscInt)(messageLength/numMessages) : 0;
15297929ea7SJunchao Zhang   /* when nproc > 2, tot_msg = 2*nproc*niter, tot_len = tot_msg*sizeof(PetscScalar)*bs */
153*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"VecScatter(bs=%" PetscInt_FMT ") has sent out %" PetscInt_FMT " messages, total %" PetscInt_FMT " bytes, with average length %" PetscInt_FMT " bytes\n",bs,tot_msg,tot_len,avg_len));
15497929ea7SJunchao Zhang #endif
15597929ea7SJunchao Zhang 
156*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Program finished\n"));
157*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isx));
158*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isy));
159*9566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&ctx));
160*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
161*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
162*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
163b122ec5aSJacob Faibussowitsch   return 0;
16497929ea7SJunchao Zhang }
16597929ea7SJunchao Zhang 
16697929ea7SJunchao Zhang /*TEST
16797929ea7SJunchao Zhang 
16897929ea7SJunchao Zhang    test:
16997929ea7SJunchao Zhang       nsize: 4
17097929ea7SJunchao Zhang       args:
171dfd57a17SPierre Jolivet       requires: double defined(PETSC_USE_LOG)
17297929ea7SJunchao Zhang 
17397929ea7SJunchao Zhang TEST*/
174