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*327415f7SBarry Smith PetscFunctionBeginUser; 279566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 289566063dSJacob Faibussowitsch PetscCall(PetscLogDefaultBegin()); 299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&nproc)); 309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 3197929ea7SJunchao Zhang 329566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Scatter(bs=1)", &stage1)); 339566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("VecScatter(bs=1)", PETSC_OBJECT_CLASSID, &event1)); 349566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Scatter(bs=4)", &stage2)); 359566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("VecScatter(bs=4)", PETSC_OBJECT_CLASSID, &event2)); 3697929ea7SJunchao Zhang 3797929ea7SJunchao Zhang /* Create a parallel vector x and a sequential vector y */ 389566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&x)); 399566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x,n,PETSC_DECIDE)); 409566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 419566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x,&low,&high)); 429566063dSJacob Faibussowitsch PetscCall(VecGetSize(x,&N)); 439566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&y)); 4497929ea7SJunchao Zhang 4597929ea7SJunchao Zhang /*======================================= 4697929ea7SJunchao Zhang test VecScatter with bs = 1 4797929ea7SJunchao Zhang ======================================*/ 4897929ea7SJunchao Zhang 4997929ea7SJunchao Zhang /* the code works as if we are going to do 3-point stencil computations on a 1D domain x, 5097929ea7SJunchao Zhang which has periodic boundary conditions but the two ghosts are scatterred to beginning of y. 5197929ea7SJunchao Zhang */ 5297929ea7SJunchao Zhang bs = 1; 5397929ea7SJunchao Zhang ix[0] = rank? low-1 : N-1; /* ix[] contains global indices of the two ghost points */ 5497929ea7SJunchao Zhang ix[1] = (rank != nproc-1)? high : 0; 5597929ea7SJunchao Zhang iy[0] = 0; 5697929ea7SJunchao Zhang iy[1] = 1; 5797929ea7SJunchao Zhang 589566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,2,ix,PETSC_COPY_VALUES,&isx)); 599566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,2,iy,PETSC_COPY_VALUES,&isy)); 609566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x,isx,y,isy,&ctx)); 619566063dSJacob Faibussowitsch PetscCall(VecScatterSetUp(ctx)); 6297929ea7SJunchao Zhang 639566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage1)); 649566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(event1,0,0,0,0)); 6597929ea7SJunchao Zhang errors = 0; 6697929ea7SJunchao Zhang for (i=0; i<niter; i++) { 6797929ea7SJunchao Zhang /* set x = 0+i, 1+i, 2+i, ..., N-1+i */ 689566063dSJacob Faibussowitsch PetscCall(VecGetArray(x,&xval)); 6997929ea7SJunchao Zhang for (j=0; j<n; j++) xval[j] = (PetscScalar)(low+j+i); 709566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x,&xval)); 7197929ea7SJunchao Zhang /* scatter the ghosts to y */ 729566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD)); 739566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD)); 7497929ea7SJunchao Zhang /* check if y has correct values */ 759566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(y,&yval)); 7697929ea7SJunchao Zhang if ((PetscInt)PetscRealPart(yval[0]) != ix[0]+i) errors++; 7797929ea7SJunchao Zhang if ((PetscInt)PetscRealPart(yval[1]) != ix[1]+i) errors++; 789566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(y,&yval)); 7997929ea7SJunchao Zhang } 809566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(event1,0,0,0,0)); 819566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 8297929ea7SJunchao Zhang 8397929ea7SJunchao Zhang /* check if we found wrong values on any processors */ 849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&errors,&tot_errors,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD)); 859566063dSJacob Faibussowitsch if (tot_errors) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error: wrong values were scatterred in vecscatter with bs = %" PetscInt_FMT "\n",bs)); 8697929ea7SJunchao Zhang 8797929ea7SJunchao Zhang /* print out event log of VecScatter(bs=1) */ 8897929ea7SJunchao Zhang #if defined(PETSC_USE_LOG) 899566063dSJacob Faibussowitsch PetscCall(PetscLogEventGetPerfInfo(stage1,event1,&eventInfo)); 909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&eventInfo.numMessages, &numMessages, 1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD)); 919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&eventInfo.messageLength,&messageLength,1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD)); 9297929ea7SJunchao Zhang tot_msg = (PetscInt)numMessages*0.5; /* two MPI calls (Send & Recv) per message */ 9397929ea7SJunchao Zhang tot_len = (PetscInt)messageLength*0.5; 9497929ea7SJunchao Zhang avg_len = tot_msg? (PetscInt)(messageLength/numMessages) : 0; 9597929ea7SJunchao Zhang /* when nproc > 2, tot_msg = 2*nproc*niter, tot_len = tot_msg*sizeof(PetscScalar)*bs */ 969566063dSJacob 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)); 9797929ea7SJunchao Zhang #endif 9897929ea7SJunchao Zhang 999566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isx)); 1009566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isy)); 1019566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ctx)); 10297929ea7SJunchao Zhang 10397929ea7SJunchao Zhang /*======================================= 10497929ea7SJunchao Zhang test VecScatter with bs = 4 10597929ea7SJunchao Zhang ======================================*/ 10697929ea7SJunchao Zhang 10797929ea7SJunchao Zhang /* similar to the 3-point stencil above, except that this time a ghost is a block */ 10897929ea7SJunchao Zhang bs = 4; /* n needs to be a multiple of bs to make the following code work */ 10997929ea7SJunchao Zhang ix[0] = rank? low/bs-1 : N/bs-1; /* ix[] contains global indices of the two ghost blocks */ 11097929ea7SJunchao Zhang ix[1] = (rank != nproc-1)? high/bs : 0; 11197929ea7SJunchao Zhang iy[0] = 0; 11297929ea7SJunchao Zhang iy[1] = 1; 11397929ea7SJunchao Zhang 1149566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,2,ix,PETSC_COPY_VALUES,&isx)); 1159566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,2,iy,PETSC_COPY_VALUES,&isy)); 11697929ea7SJunchao Zhang 1179566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x,isx,y,isy,&ctx)); 11897929ea7SJunchao Zhang /* Call SetUp explicitly, otherwise messages in implicit SetUp will be counted in events below */ 1199566063dSJacob Faibussowitsch PetscCall(VecScatterSetUp(ctx)); 12097929ea7SJunchao Zhang 1219566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage2)); 1229566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(event2,0,0,0,0)); 12397929ea7SJunchao Zhang errors = 0; 12497929ea7SJunchao Zhang for (i=0; i<niter; i++) { 12597929ea7SJunchao Zhang /* set x = 0+i, 1+i, 2+i, ..., N-1+i */ 1269566063dSJacob Faibussowitsch PetscCall(VecGetArray(x,&xval)); 12797929ea7SJunchao Zhang for (j=0; j<n; j++) xval[j] = (PetscScalar)(low+j+i); 1289566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x,&xval)); 12997929ea7SJunchao Zhang /* scatter the ghost blocks to y */ 1309566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD)); 1319566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD)); 13297929ea7SJunchao Zhang /* check if y has correct values */ 1339566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(y,&yval)); 13497929ea7SJunchao Zhang if ((PetscInt)PetscRealPart(yval[0]) != ix[0]*bs+i) errors++; 13597929ea7SJunchao Zhang if ((PetscInt)PetscRealPart(yval[bs]) != ix[1]*bs+i) errors++; 1369566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(y,&yval)); 13797929ea7SJunchao Zhang } 1389566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(event2,0,0,0,0)); 1399566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 14097929ea7SJunchao Zhang 14197929ea7SJunchao Zhang /* check if we found wrong values on any processors */ 1429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&errors,&tot_errors,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD)); 1439566063dSJacob Faibussowitsch if (tot_errors) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error: wrong values were scatterred in vecscatter with bs = %" PetscInt_FMT "\n",bs)); 14497929ea7SJunchao Zhang 14597929ea7SJunchao Zhang /* print out event log of VecScatter(bs=4) */ 14697929ea7SJunchao Zhang #if defined(PETSC_USE_LOG) 1479566063dSJacob Faibussowitsch PetscCall(PetscLogEventGetPerfInfo(stage2,event2,&eventInfo)); 1489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&eventInfo.numMessages, &numMessages, 1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD)); 1499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&eventInfo.messageLength,&messageLength,1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD)); 15097929ea7SJunchao Zhang tot_msg = (PetscInt)numMessages*0.5; /* two MPI calls (Send & Recv) per message */ 15197929ea7SJunchao Zhang tot_len = (PetscInt)messageLength*0.5; 15297929ea7SJunchao Zhang avg_len = tot_msg? (PetscInt)(messageLength/numMessages) : 0; 15397929ea7SJunchao Zhang /* when nproc > 2, tot_msg = 2*nproc*niter, tot_len = tot_msg*sizeof(PetscScalar)*bs */ 1549566063dSJacob 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)); 15597929ea7SJunchao Zhang #endif 15697929ea7SJunchao Zhang 1579566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Program finished\n")); 1589566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isx)); 1599566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isy)); 1609566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ctx)); 1619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1639566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 164b122ec5aSJacob Faibussowitsch return 0; 16597929ea7SJunchao Zhang } 16697929ea7SJunchao Zhang 16797929ea7SJunchao Zhang /*TEST 16897929ea7SJunchao Zhang 16997929ea7SJunchao Zhang test: 17097929ea7SJunchao Zhang nsize: 4 17197929ea7SJunchao Zhang args: 172dfd57a17SPierre Jolivet requires: double defined(PETSC_USE_LOG) 17397929ea7SJunchao Zhang 17497929ea7SJunchao Zhang TEST*/ 175