xref: /petsc/src/vec/is/sf/tests/ex14.c (revision 2611ad710242b3c70d66651ef7e40f9450d305e2)
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 
6d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
7d71ae5a4SJacob Faibussowitsch {
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   PetscLogStage      stage1, stage2;
1897929ea7SJunchao Zhang   PetscLogEvent      event1, event2;
1997929ea7SJunchao Zhang 
2097929ea7SJunchao Zhang   PetscFunctionBegin;
21327415f7SBarry Smith   PetscFunctionBeginUser;
229566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
239566063dSJacob Faibussowitsch   PetscCall(PetscLogDefaultBegin());
249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &nproc));
259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
2697929ea7SJunchao Zhang 
279566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("Scatter(bs=1)", &stage1));
289566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("VecScatter(bs=1)", PETSC_OBJECT_CLASSID, &event1));
299566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("Scatter(bs=4)", &stage2));
309566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("VecScatter(bs=4)", PETSC_OBJECT_CLASSID, &event2));
3197929ea7SJunchao Zhang 
3297929ea7SJunchao Zhang   /* Create a parallel vector x and a sequential vector y */
339566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
349566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, n, PETSC_DECIDE));
359566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
369566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x, &low, &high));
379566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &N));
389566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &y));
3997929ea7SJunchao Zhang 
4097929ea7SJunchao Zhang   /*=======================================
4197929ea7SJunchao Zhang      test VecScatter with bs = 1
4297929ea7SJunchao Zhang     ======================================*/
4397929ea7SJunchao Zhang 
4497929ea7SJunchao Zhang   /* the code works as if we are going to do 3-point stencil computations on a 1D domain x,
4597929ea7SJunchao Zhang      which has periodic boundary conditions but the two ghosts are scatterred to beginning of y.
4697929ea7SJunchao Zhang    */
4797929ea7SJunchao Zhang   bs    = 1;
4897929ea7SJunchao Zhang   ix[0] = rank ? low - 1 : N - 1; /* ix[] contains global indices of the two ghost points */
4997929ea7SJunchao Zhang   ix[1] = (rank != nproc - 1) ? high : 0;
5097929ea7SJunchao Zhang   iy[0] = 0;
5197929ea7SJunchao Zhang   iy[1] = 1;
5297929ea7SJunchao Zhang 
539566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 2, ix, PETSC_COPY_VALUES, &isx));
549566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 2, iy, PETSC_COPY_VALUES, &isy));
559566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x, isx, y, isy, &ctx));
569566063dSJacob Faibussowitsch   PetscCall(VecScatterSetUp(ctx));
5797929ea7SJunchao Zhang 
589566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage1));
599566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(event1, 0, 0, 0, 0));
6097929ea7SJunchao Zhang   errors = 0;
6197929ea7SJunchao Zhang   for (i = 0; i < niter; i++) {
6297929ea7SJunchao Zhang     /* set x = 0+i, 1+i, 2+i, ..., N-1+i */
639566063dSJacob Faibussowitsch     PetscCall(VecGetArray(x, &xval));
6497929ea7SJunchao Zhang     for (j = 0; j < n; j++) xval[j] = (PetscScalar)(low + j + i);
659566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(x, &xval));
6697929ea7SJunchao Zhang     /* scatter the ghosts to y */
679566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ctx, x, y, INSERT_VALUES, SCATTER_FORWARD));
689566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ctx, x, y, INSERT_VALUES, SCATTER_FORWARD));
6997929ea7SJunchao Zhang     /* check if y has correct values */
709566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(y, &yval));
7197929ea7SJunchao Zhang     if ((PetscInt)PetscRealPart(yval[0]) != ix[0] + i) errors++;
7297929ea7SJunchao Zhang     if ((PetscInt)PetscRealPart(yval[1]) != ix[1] + i) errors++;
739566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(y, &yval));
7497929ea7SJunchao Zhang   }
759566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(event1, 0, 0, 0, 0));
769566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
7797929ea7SJunchao Zhang 
7897929ea7SJunchao Zhang   /* check if we found wrong values on any processors */
79712fec58SPierre Jolivet   PetscCall(MPIU_Allreduce(&errors, &tot_errors, 1, MPIU_INT, MPI_SUM, PETSC_COMM_WORLD));
809566063dSJacob Faibussowitsch   if (tot_errors) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: wrong values were scatterred in vecscatter with bs = %" PetscInt_FMT "\n", bs));
8197929ea7SJunchao Zhang 
8297929ea7SJunchao Zhang   /* print out event log of VecScatter(bs=1) */
83*2611ad71SToby Isaac   if (PetscDefined(USE_LOG)) {
84*2611ad71SToby Isaac     PetscLogDouble     numMessages, messageLength;
85*2611ad71SToby Isaac     PetscEventPerfInfo eventInfo;
86*2611ad71SToby Isaac     PetscInt           tot_msg, tot_len, avg_len;
87*2611ad71SToby Isaac 
889566063dSJacob Faibussowitsch     PetscCall(PetscLogEventGetPerfInfo(stage1, event1, &eventInfo));
89712fec58SPierre Jolivet     PetscCall(MPIU_Allreduce(&eventInfo.numMessages, &numMessages, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, PETSC_COMM_WORLD));
90712fec58SPierre Jolivet     PetscCall(MPIU_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 */
959566063dSJacob 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));
96*2611ad71SToby Isaac   }
9797929ea7SJunchao Zhang 
989566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isx));
999566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isy));
1009566063dSJacob 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 
1139566063dSJacob Faibussowitsch   PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 2, ix, PETSC_COPY_VALUES, &isx));
1149566063dSJacob Faibussowitsch   PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 2, iy, PETSC_COPY_VALUES, &isy));
11597929ea7SJunchao Zhang 
1169566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x, isx, y, isy, &ctx));
11797929ea7SJunchao Zhang   /* Call SetUp explicitly, otherwise messages in implicit SetUp will be counted in events below */
1189566063dSJacob Faibussowitsch   PetscCall(VecScatterSetUp(ctx));
11997929ea7SJunchao Zhang 
1209566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage2));
1219566063dSJacob 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 */
1259566063dSJacob Faibussowitsch     PetscCall(VecGetArray(x, &xval));
12697929ea7SJunchao Zhang     for (j = 0; j < n; j++) xval[j] = (PetscScalar)(low + j + i);
1279566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(x, &xval));
12897929ea7SJunchao Zhang     /* scatter the ghost blocks to y */
1299566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ctx, x, y, INSERT_VALUES, SCATTER_FORWARD));
1309566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(ctx, x, y, INSERT_VALUES, SCATTER_FORWARD));
13197929ea7SJunchao Zhang     /* check if y has correct values */
1329566063dSJacob 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++;
1359566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(y, &yval));
13697929ea7SJunchao Zhang   }
1379566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(event2, 0, 0, 0, 0));
1389566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
13997929ea7SJunchao Zhang 
14097929ea7SJunchao Zhang   /* check if we found wrong values on any processors */
141712fec58SPierre Jolivet   PetscCall(MPIU_Allreduce(&errors, &tot_errors, 1, MPIU_INT, MPI_SUM, PETSC_COMM_WORLD));
1429566063dSJacob 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) */
145*2611ad71SToby Isaac   if (PetscDefined(USE_LOG)) {
146*2611ad71SToby Isaac     PetscLogDouble     numMessages, messageLength;
147*2611ad71SToby Isaac     PetscEventPerfInfo eventInfo;
148*2611ad71SToby Isaac     PetscInt           tot_msg, tot_len, avg_len;
149*2611ad71SToby Isaac 
1509566063dSJacob Faibussowitsch     PetscCall(PetscLogEventGetPerfInfo(stage2, event2, &eventInfo));
151712fec58SPierre Jolivet     PetscCall(MPIU_Allreduce(&eventInfo.numMessages, &numMessages, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, PETSC_COMM_WORLD));
152712fec58SPierre Jolivet     PetscCall(MPIU_Allreduce(&eventInfo.messageLength, &messageLength, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, PETSC_COMM_WORLD));
15397929ea7SJunchao Zhang     tot_msg = (PetscInt)numMessages * 0.5; /* two MPI calls (Send & Recv) per message */
15497929ea7SJunchao Zhang     tot_len = (PetscInt)messageLength * 0.5;
15597929ea7SJunchao Zhang     avg_len = tot_msg ? (PetscInt)(messageLength / numMessages) : 0;
15697929ea7SJunchao Zhang     /* when nproc > 2, tot_msg = 2*nproc*niter, tot_len = tot_msg*sizeof(PetscScalar)*bs */
1579566063dSJacob 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));
158*2611ad71SToby Isaac   }
15997929ea7SJunchao Zhang 
1609566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Program finished\n"));
1619566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isx));
1629566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isy));
1639566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&ctx));
1649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
1669566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
167b122ec5aSJacob Faibussowitsch   return 0;
16897929ea7SJunchao Zhang }
16997929ea7SJunchao Zhang 
17097929ea7SJunchao Zhang /*TEST
17197929ea7SJunchao Zhang 
17297929ea7SJunchao Zhang    test:
17397929ea7SJunchao Zhang       nsize: 4
17497929ea7SJunchao Zhang       args:
175dfd57a17SPierre Jolivet       requires: double defined(PETSC_USE_LOG)
17697929ea7SJunchao Zhang 
17797929ea7SJunchao Zhang TEST*/
178