197929ea7SJunchao Zhang static char help[] = "Test event log of VecScatter with various block sizes\n\n"; 297929ea7SJunchao Zhang 397929ea7SJunchao Zhang #include <petscvec.h> 497929ea7SJunchao Zhang 5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 6d71ae5a4SJacob Faibussowitsch { 797929ea7SJunchao Zhang PetscInt i, j, low, high, n = 256, N, errors, tot_errors; 897929ea7SJunchao Zhang PetscInt bs = 1, ix[2], iy[2]; 997929ea7SJunchao Zhang PetscMPIInt nproc, rank; 1097929ea7SJunchao Zhang PetscScalar *xval; 1197929ea7SJunchao Zhang const PetscScalar *yval; 1297929ea7SJunchao Zhang Vec x, y; 1397929ea7SJunchao Zhang IS isx, isy; 1497929ea7SJunchao Zhang VecScatter ctx; 1597929ea7SJunchao Zhang const PetscInt niter = 10; 1697929ea7SJunchao Zhang PetscLogStage stage1, stage2; 1797929ea7SJunchao Zhang PetscLogEvent event1, event2; 1897929ea7SJunchao Zhang 1997929ea7SJunchao Zhang PetscFunctionBegin; 20327415f7SBarry Smith PetscFunctionBeginUser; 21*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 229566063dSJacob Faibussowitsch PetscCall(PetscLogDefaultBegin()); 239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &nproc)); 249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 2597929ea7SJunchao Zhang 269566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Scatter(bs=1)", &stage1)); 279566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("VecScatter(bs=1)", PETSC_OBJECT_CLASSID, &event1)); 289566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Scatter(bs=4)", &stage2)); 299566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("VecScatter(bs=4)", PETSC_OBJECT_CLASSID, &event2)); 3097929ea7SJunchao Zhang 3197929ea7SJunchao Zhang /* Create a parallel vector x and a sequential vector y */ 329566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 339566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, n, PETSC_DECIDE)); 349566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 359566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x, &low, &high)); 369566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &N)); 379566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &y)); 3897929ea7SJunchao Zhang 3997929ea7SJunchao Zhang /*======================================= 4097929ea7SJunchao Zhang test VecScatter with bs = 1 4197929ea7SJunchao Zhang ======================================*/ 4297929ea7SJunchao Zhang 4397929ea7SJunchao Zhang /* the code works as if we are going to do 3-point stencil computations on a 1D domain x, 4497929ea7SJunchao Zhang which has periodic boundary conditions but the two ghosts are scatterred to beginning of y. 4597929ea7SJunchao Zhang */ 4697929ea7SJunchao Zhang bs = 1; 4797929ea7SJunchao Zhang ix[0] = rank ? low - 1 : N - 1; /* ix[] contains global indices of the two ghost points */ 4897929ea7SJunchao Zhang ix[1] = (rank != nproc - 1) ? high : 0; 4997929ea7SJunchao Zhang iy[0] = 0; 5097929ea7SJunchao Zhang iy[1] = 1; 5197929ea7SJunchao Zhang 529566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 2, ix, PETSC_COPY_VALUES, &isx)); 539566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 2, iy, PETSC_COPY_VALUES, &isy)); 549566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x, isx, y, isy, &ctx)); 559566063dSJacob Faibussowitsch PetscCall(VecScatterSetUp(ctx)); 5697929ea7SJunchao Zhang 579566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage1)); 589566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(event1, 0, 0, 0, 0)); 5997929ea7SJunchao Zhang errors = 0; 6097929ea7SJunchao Zhang for (i = 0; i < niter; i++) { 6197929ea7SJunchao Zhang /* set x = 0+i, 1+i, 2+i, ..., N-1+i */ 629566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &xval)); 6397929ea7SJunchao Zhang for (j = 0; j < n; j++) xval[j] = (PetscScalar)(low + j + i); 649566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &xval)); 6597929ea7SJunchao Zhang /* scatter the ghosts to y */ 669566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx, x, y, INSERT_VALUES, SCATTER_FORWARD)); 679566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx, x, y, INSERT_VALUES, SCATTER_FORWARD)); 6897929ea7SJunchao Zhang /* check if y has correct values */ 699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(y, &yval)); 7097929ea7SJunchao Zhang if ((PetscInt)PetscRealPart(yval[0]) != ix[0] + i) errors++; 7197929ea7SJunchao Zhang if ((PetscInt)PetscRealPart(yval[1]) != ix[1] + i) errors++; 729566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(y, &yval)); 7397929ea7SJunchao Zhang } 749566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(event1, 0, 0, 0, 0)); 759566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 7697929ea7SJunchao Zhang 7797929ea7SJunchao Zhang /* check if we found wrong values on any processors */ 78462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&errors, &tot_errors, 1, MPIU_INT, MPI_SUM, PETSC_COMM_WORLD)); 799566063dSJacob Faibussowitsch if (tot_errors) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: wrong values were scatterred in vecscatter with bs = %" PetscInt_FMT "\n", bs)); 8097929ea7SJunchao Zhang 8197929ea7SJunchao Zhang /* print out event log of VecScatter(bs=1) */ 822611ad71SToby Isaac if (PetscDefined(USE_LOG)) { 832611ad71SToby Isaac PetscLogDouble numMessages, messageLength; 842611ad71SToby Isaac PetscEventPerfInfo eventInfo; 852611ad71SToby Isaac PetscInt tot_msg, tot_len, avg_len; 862611ad71SToby Isaac 879566063dSJacob Faibussowitsch PetscCall(PetscLogEventGetPerfInfo(stage1, event1, &eventInfo)); 88462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&eventInfo.numMessages, &numMessages, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, PETSC_COMM_WORLD)); 89462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&eventInfo.messageLength, &messageLength, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, PETSC_COMM_WORLD)); 9097929ea7SJunchao Zhang tot_msg = (PetscInt)numMessages * 0.5; /* two MPI calls (Send & Recv) per message */ 9197929ea7SJunchao Zhang tot_len = (PetscInt)messageLength * 0.5; 9297929ea7SJunchao Zhang avg_len = tot_msg ? (PetscInt)(messageLength / numMessages) : 0; 9397929ea7SJunchao Zhang /* when nproc > 2, tot_msg = 2*nproc*niter, tot_len = tot_msg*sizeof(PetscScalar)*bs */ 949566063dSJacob 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)); 952611ad71SToby Isaac } 9697929ea7SJunchao Zhang 979566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isx)); 989566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isy)); 999566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ctx)); 10097929ea7SJunchao Zhang 10197929ea7SJunchao Zhang /*======================================= 10297929ea7SJunchao Zhang test VecScatter with bs = 4 10397929ea7SJunchao Zhang ======================================*/ 10497929ea7SJunchao Zhang 10597929ea7SJunchao Zhang /* similar to the 3-point stencil above, except that this time a ghost is a block */ 10697929ea7SJunchao Zhang bs = 4; /* n needs to be a multiple of bs to make the following code work */ 10797929ea7SJunchao Zhang ix[0] = rank ? low / bs - 1 : N / bs - 1; /* ix[] contains global indices of the two ghost blocks */ 10897929ea7SJunchao Zhang ix[1] = (rank != nproc - 1) ? high / bs : 0; 10997929ea7SJunchao Zhang iy[0] = 0; 11097929ea7SJunchao Zhang iy[1] = 1; 11197929ea7SJunchao Zhang 1129566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 2, ix, PETSC_COPY_VALUES, &isx)); 1139566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 2, iy, PETSC_COPY_VALUES, &isy)); 11497929ea7SJunchao Zhang 1159566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x, isx, y, isy, &ctx)); 11697929ea7SJunchao Zhang /* Call SetUp explicitly, otherwise messages in implicit SetUp will be counted in events below */ 1179566063dSJacob Faibussowitsch PetscCall(VecScatterSetUp(ctx)); 11897929ea7SJunchao Zhang 1199566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage2)); 1209566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(event2, 0, 0, 0, 0)); 12197929ea7SJunchao Zhang errors = 0; 12297929ea7SJunchao Zhang for (i = 0; i < niter; i++) { 12397929ea7SJunchao Zhang /* set x = 0+i, 1+i, 2+i, ..., N-1+i */ 1249566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &xval)); 12597929ea7SJunchao Zhang for (j = 0; j < n; j++) xval[j] = (PetscScalar)(low + j + i); 1269566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &xval)); 12797929ea7SJunchao Zhang /* scatter the ghost blocks to y */ 1289566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx, x, y, INSERT_VALUES, SCATTER_FORWARD)); 1299566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx, x, y, INSERT_VALUES, SCATTER_FORWARD)); 13097929ea7SJunchao Zhang /* check if y has correct values */ 1319566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(y, &yval)); 13297929ea7SJunchao Zhang if ((PetscInt)PetscRealPart(yval[0]) != ix[0] * bs + i) errors++; 13397929ea7SJunchao Zhang if ((PetscInt)PetscRealPart(yval[bs]) != ix[1] * bs + i) errors++; 1349566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(y, &yval)); 13597929ea7SJunchao Zhang } 1369566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(event2, 0, 0, 0, 0)); 1379566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 13897929ea7SJunchao Zhang 13997929ea7SJunchao Zhang /* check if we found wrong values on any processors */ 140462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&errors, &tot_errors, 1, MPIU_INT, MPI_SUM, PETSC_COMM_WORLD)); 1419566063dSJacob Faibussowitsch if (tot_errors) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: wrong values were scatterred in vecscatter with bs = %" PetscInt_FMT "\n", bs)); 14297929ea7SJunchao Zhang 14397929ea7SJunchao Zhang /* print out event log of VecScatter(bs=4) */ 1442611ad71SToby Isaac if (PetscDefined(USE_LOG)) { 1452611ad71SToby Isaac PetscLogDouble numMessages, messageLength; 1462611ad71SToby Isaac PetscEventPerfInfo eventInfo; 1472611ad71SToby Isaac PetscInt tot_msg, tot_len, avg_len; 1482611ad71SToby Isaac 1499566063dSJacob Faibussowitsch PetscCall(PetscLogEventGetPerfInfo(stage2, event2, &eventInfo)); 150462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&eventInfo.numMessages, &numMessages, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, PETSC_COMM_WORLD)); 151462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&eventInfo.messageLength, &messageLength, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, PETSC_COMM_WORLD)); 15297929ea7SJunchao Zhang tot_msg = (PetscInt)numMessages * 0.5; /* two MPI calls (Send & Recv) per message */ 15397929ea7SJunchao Zhang tot_len = (PetscInt)messageLength * 0.5; 15497929ea7SJunchao Zhang avg_len = tot_msg ? (PetscInt)(messageLength / numMessages) : 0; 15597929ea7SJunchao Zhang /* when nproc > 2, tot_msg = 2*nproc*niter, tot_len = tot_msg*sizeof(PetscScalar)*bs */ 1569566063dSJacob 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)); 1572611ad71SToby Isaac } 15897929ea7SJunchao Zhang 1599566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Program finished\n")); 1609566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isx)); 1619566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isy)); 1629566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ctx)); 1639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1659566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 166b122ec5aSJacob Faibussowitsch return 0; 16797929ea7SJunchao Zhang } 16897929ea7SJunchao Zhang 16997929ea7SJunchao Zhang /*TEST 17097929ea7SJunchao Zhang 17197929ea7SJunchao Zhang test: 17297929ea7SJunchao Zhang nsize: 4 17397929ea7SJunchao Zhang args: 174dfd57a17SPierre Jolivet requires: double defined(PETSC_USE_LOG) 17597929ea7SJunchao Zhang 17697929ea7SJunchao Zhang TEST*/ 177