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*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogDefaultBegin()); 285f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&nproc)); 295f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 3097929ea7SJunchao Zhang 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister("Scatter(bs=1)", &stage1)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventRegister("VecScatter(bs=1)", PETSC_OBJECT_CLASSID, &event1)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister("Scatter(bs=4)", &stage2)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventRegister("VecScatter(bs=4)", PETSC_OBJECT_CLASSID, &event2)); 3597929ea7SJunchao Zhang 3697929ea7SJunchao Zhang /* Create a parallel vector x and a sequential vector y */ 375f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(x,n,PETSC_DECIDE)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(x)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(x,&low,&high)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(x,&N)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(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 575f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,2,ix,PETSC_COPY_VALUES,&isx)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,2,iy,PETSC_COPY_VALUES,&isy)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(x,isx,y,isy,&ctx)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterSetUp(ctx)); 6197929ea7SJunchao Zhang 625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(stage1)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 675f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(x,&xval)); 6897929ea7SJunchao Zhang for (j=0; j<n; j++) xval[j] = (PetscScalar)(low+j+i); 695f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(x,&xval)); 7097929ea7SJunchao Zhang /* scatter the ghosts to y */ 715f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD)); 7397929ea7SJunchao Zhang /* check if y has correct values */ 745f80ce2aSJacob Faibussowitsch CHKERRQ(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++; 775f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(y,&yval)); 7897929ea7SJunchao Zhang } 795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventEnd(event1,0,0,0,0)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 8197929ea7SJunchao Zhang 8297929ea7SJunchao Zhang /* check if we found wrong values on any processors */ 835f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(&errors,&tot_errors,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD)); 845f80ce2aSJacob Faibussowitsch if (tot_errors) CHKERRQ(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) 885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventGetPerfInfo(stage1,event1,&eventInfo)); 895f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(&eventInfo.numMessages, &numMessages, 1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD)); 905f80ce2aSJacob Faibussowitsch CHKERRMPI(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 */ 955f80ce2aSJacob Faibussowitsch CHKERRQ(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 985f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isx)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isy)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(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 1135f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,2,ix,PETSC_COPY_VALUES,&isx)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,2,iy,PETSC_COPY_VALUES,&isy)); 11597929ea7SJunchao Zhang 1165f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(x,isx,y,isy,&ctx)); 11797929ea7SJunchao Zhang /* Call SetUp explicitly, otherwise messages in implicit SetUp will be counted in events below */ 1185f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterSetUp(ctx)); 11997929ea7SJunchao Zhang 1205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(stage2)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 1255f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(x,&xval)); 12697929ea7SJunchao Zhang for (j=0; j<n; j++) xval[j] = (PetscScalar)(low+j+i); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(x,&xval)); 12897929ea7SJunchao Zhang /* scatter the ghost blocks to y */ 1295f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD)); 13197929ea7SJunchao Zhang /* check if y has correct values */ 1325f80ce2aSJacob Faibussowitsch CHKERRQ(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++; 1355f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(y,&yval)); 13697929ea7SJunchao Zhang } 1375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventEnd(event2,0,0,0,0)); 1385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 13997929ea7SJunchao Zhang 14097929ea7SJunchao Zhang /* check if we found wrong values on any processors */ 1415f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(&errors,&tot_errors,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD)); 1425f80ce2aSJacob Faibussowitsch if (tot_errors) CHKERRQ(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) 1465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventGetPerfInfo(stage2,event2,&eventInfo)); 1475f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(&eventInfo.numMessages, &numMessages, 1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD)); 1485f80ce2aSJacob Faibussowitsch CHKERRMPI(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 */ 1535f80ce2aSJacob Faibussowitsch CHKERRQ(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 1565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Program finished\n")); 1575f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isx)); 1585f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isy)); 1595f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&ctx)); 1605f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 1615f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 162*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 163*b122ec5aSJacob 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