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 PetscErrorCode ierr; 997929ea7SJunchao Zhang PetscInt i,j,low,high,n=256,N,errors,tot_errors; 1097929ea7SJunchao Zhang PetscInt bs=1,ix[2],iy[2]; 1197929ea7SJunchao Zhang PetscMPIInt nproc,rank; 1297929ea7SJunchao Zhang PetscScalar *xval; 1397929ea7SJunchao Zhang const PetscScalar *yval; 1497929ea7SJunchao Zhang Vec x,y; 1597929ea7SJunchao Zhang IS isx,isy; 1697929ea7SJunchao Zhang VecScatter ctx; 1797929ea7SJunchao Zhang const PetscInt niter = 10; 1897929ea7SJunchao Zhang #if defined(PETSC_USE_LOG) 1997929ea7SJunchao Zhang PetscLogStage stage1,stage2; 2097929ea7SJunchao Zhang PetscLogEvent event1,event2; 2197929ea7SJunchao Zhang PetscLogDouble numMessages,messageLength; 2297929ea7SJunchao Zhang PetscEventPerfInfo eventInfo; 2397929ea7SJunchao Zhang PetscInt tot_msg,tot_len,avg_len; 2497929ea7SJunchao Zhang #endif 2597929ea7SJunchao Zhang 2697929ea7SJunchao Zhang PetscFunctionBegin; 2797929ea7SJunchao Zhang ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 2897929ea7SJunchao Zhang ierr = PetscLogDefaultBegin();CHKERRQ(ierr); 29ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&nproc);CHKERRMPI(ierr); 30ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 3197929ea7SJunchao Zhang 3297929ea7SJunchao Zhang ierr = PetscLogStageRegister("Scatter(bs=1)", &stage1);CHKERRQ(ierr); 3397929ea7SJunchao Zhang ierr = PetscLogEventRegister("VecScatter(bs=1)", PETSC_OBJECT_CLASSID, &event1);CHKERRQ(ierr); 3497929ea7SJunchao Zhang ierr = PetscLogStageRegister("Scatter(bs=4)", &stage2);CHKERRQ(ierr); 3597929ea7SJunchao Zhang ierr = PetscLogEventRegister("VecScatter(bs=4)", PETSC_OBJECT_CLASSID, &event2);CHKERRQ(ierr); 3697929ea7SJunchao Zhang 3797929ea7SJunchao Zhang /* Create a parallel vector x and a sequential vector y */ 3897929ea7SJunchao Zhang ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 3997929ea7SJunchao Zhang ierr = VecSetSizes(x,n,PETSC_DECIDE);CHKERRQ(ierr); 4097929ea7SJunchao Zhang ierr = VecSetFromOptions(x);CHKERRQ(ierr); 4197929ea7SJunchao Zhang ierr = VecGetOwnershipRange(x,&low,&high);CHKERRQ(ierr); 4297929ea7SJunchao Zhang ierr = VecGetSize(x,&N);CHKERRQ(ierr); 4397929ea7SJunchao Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,n,&y);CHKERRQ(ierr); 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 5897929ea7SJunchao Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,2,ix,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr); 5997929ea7SJunchao Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,2,iy,PETSC_COPY_VALUES,&isy);CHKERRQ(ierr); 6097929ea7SJunchao Zhang ierr = VecScatterCreate(x,isx,y,isy,&ctx);CHKERRQ(ierr); 6197929ea7SJunchao Zhang ierr = VecScatterSetUp(ctx);CHKERRQ(ierr); 6297929ea7SJunchao Zhang 6397929ea7SJunchao Zhang ierr = PetscLogStagePush(stage1);CHKERRQ(ierr); 6497929ea7SJunchao Zhang ierr = PetscLogEventBegin(event1,0,0,0,0);CHKERRQ(ierr); 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 */ 6897929ea7SJunchao Zhang ierr = VecGetArray(x,&xval);CHKERRQ(ierr); 6997929ea7SJunchao Zhang for (j=0; j<n; j++) xval[j] = (PetscScalar)(low+j+i); 7097929ea7SJunchao Zhang ierr = VecRestoreArray(x,&xval);CHKERRQ(ierr); 7197929ea7SJunchao Zhang /* scatter the ghosts to y */ 7297929ea7SJunchao Zhang ierr = VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7397929ea7SJunchao Zhang ierr = VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7497929ea7SJunchao Zhang /* check if y has correct values */ 7597929ea7SJunchao Zhang ierr = VecGetArrayRead(y,&yval);CHKERRQ(ierr); 7697929ea7SJunchao Zhang if ((PetscInt)PetscRealPart(yval[0]) != ix[0]+i) errors++; 7797929ea7SJunchao Zhang if ((PetscInt)PetscRealPart(yval[1]) != ix[1]+i) errors++; 7897929ea7SJunchao Zhang ierr = VecRestoreArrayRead(y,&yval);CHKERRQ(ierr); 7997929ea7SJunchao Zhang } 8097929ea7SJunchao Zhang ierr = PetscLogEventEnd(event1,0,0,0,0);CHKERRQ(ierr); 8197929ea7SJunchao Zhang ierr = PetscLogStagePop();CHKERRQ(ierr); 8297929ea7SJunchao Zhang 8397929ea7SJunchao Zhang /* check if we found wrong values on any processors */ 84ffc4695bSBarry Smith ierr = MPI_Allreduce(&errors,&tot_errors,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr); 8597929ea7SJunchao Zhang if (tot_errors) { ierr = PetscPrintf(PETSC_COMM_WORLD,"Error: wrong values were scatterred in vecscatter with bs = %d\n",bs);CHKERRQ(ierr); } 8697929ea7SJunchao Zhang 8797929ea7SJunchao Zhang /* print out event log of VecScatter(bs=1) */ 8897929ea7SJunchao Zhang #if defined(PETSC_USE_LOG) 8997929ea7SJunchao Zhang ierr = PetscLogEventGetPerfInfo(stage1,event1,&eventInfo);CHKERRQ(ierr); 9055b25c41SPierre Jolivet ierr = MPI_Allreduce(&eventInfo.numMessages, &numMessages, 1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr); 9155b25c41SPierre Jolivet ierr = MPI_Allreduce(&eventInfo.messageLength,&messageLength,1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr); 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 */ 9697929ea7SJunchao Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"VecScatter(bs=%d) has sent out %d messages, total %d bytes, with average length %d bytes\n",bs,tot_msg,tot_len,avg_len);CHKERRQ(ierr); 9797929ea7SJunchao Zhang #endif 9897929ea7SJunchao Zhang 9997929ea7SJunchao Zhang ierr = ISDestroy(&isx);CHKERRQ(ierr); 10097929ea7SJunchao Zhang ierr = ISDestroy(&isy);CHKERRQ(ierr); 10197929ea7SJunchao Zhang ierr = VecScatterDestroy(&ctx);CHKERRQ(ierr); 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 11497929ea7SJunchao Zhang ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2,ix,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr); 11597929ea7SJunchao Zhang ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2,iy,PETSC_COPY_VALUES,&isy);CHKERRQ(ierr); 11697929ea7SJunchao Zhang 11797929ea7SJunchao Zhang ierr = VecScatterCreate(x,isx,y,isy,&ctx);CHKERRQ(ierr); 11897929ea7SJunchao Zhang /* Call SetUp explicitly, otherwise messages in implicit SetUp will be counted in events below */ 11997929ea7SJunchao Zhang ierr = VecScatterSetUp(ctx);CHKERRQ(ierr); 12097929ea7SJunchao Zhang 12197929ea7SJunchao Zhang ierr = PetscLogStagePush(stage2);CHKERRQ(ierr); 12297929ea7SJunchao Zhang ierr = PetscLogEventBegin(event2,0,0,0,0);CHKERRQ(ierr); 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 */ 12697929ea7SJunchao Zhang ierr = VecGetArray(x,&xval);CHKERRQ(ierr); 12797929ea7SJunchao Zhang for (j=0; j<n; j++) xval[j] = (PetscScalar)(low+j+i); 12897929ea7SJunchao Zhang ierr = VecRestoreArray(x,&xval);CHKERRQ(ierr); 12997929ea7SJunchao Zhang /* scatter the ghost blocks to y */ 13097929ea7SJunchao Zhang ierr = VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 13197929ea7SJunchao Zhang ierr = VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 13297929ea7SJunchao Zhang /* check if y has correct values */ 13397929ea7SJunchao Zhang ierr = VecGetArrayRead(y,&yval);CHKERRQ(ierr); 13497929ea7SJunchao Zhang if ((PetscInt)PetscRealPart(yval[0]) != ix[0]*bs+i) errors++; 13597929ea7SJunchao Zhang if ((PetscInt)PetscRealPart(yval[bs]) != ix[1]*bs+i) errors++; 13697929ea7SJunchao Zhang ierr = VecRestoreArrayRead(y,&yval);CHKERRQ(ierr); 13797929ea7SJunchao Zhang } 13897929ea7SJunchao Zhang ierr = PetscLogEventEnd(event2,0,0,0,0);CHKERRQ(ierr); 13997929ea7SJunchao Zhang ierr = PetscLogStagePop();CHKERRQ(ierr); 14097929ea7SJunchao Zhang 14197929ea7SJunchao Zhang /* check if we found wrong values on any processors */ 142ffc4695bSBarry Smith ierr = MPI_Allreduce(&errors,&tot_errors,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr); 14397929ea7SJunchao Zhang if (tot_errors) { ierr = PetscPrintf(PETSC_COMM_WORLD,"Error: wrong values were scatterred in vecscatter with bs = %d\n",bs);CHKERRQ(ierr); } 14497929ea7SJunchao Zhang 14597929ea7SJunchao Zhang /* print out event log of VecScatter(bs=4) */ 14697929ea7SJunchao Zhang #if defined(PETSC_USE_LOG) 14797929ea7SJunchao Zhang ierr = PetscLogEventGetPerfInfo(stage2,event2,&eventInfo);CHKERRQ(ierr); 14855b25c41SPierre Jolivet ierr = MPI_Allreduce(&eventInfo.numMessages, &numMessages, 1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr); 14955b25c41SPierre Jolivet ierr = MPI_Allreduce(&eventInfo.messageLength,&messageLength,1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr); 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 */ 15497929ea7SJunchao Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"VecScatter(bs=%d) has sent out %d messages, total %d bytes, with average length %d bytes\n",bs,tot_msg,tot_len,avg_len);CHKERRQ(ierr); 15597929ea7SJunchao Zhang #endif 15697929ea7SJunchao Zhang 15797929ea7SJunchao Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"Program finished\n");CHKERRQ(ierr); 15897929ea7SJunchao Zhang ierr = ISDestroy(&isx);CHKERRQ(ierr); 15997929ea7SJunchao Zhang ierr = ISDestroy(&isy);CHKERRQ(ierr); 16097929ea7SJunchao Zhang ierr = VecScatterDestroy(&ctx);CHKERRQ(ierr); 16197929ea7SJunchao Zhang ierr = VecDestroy(&x);CHKERRQ(ierr); 16297929ea7SJunchao Zhang ierr = VecDestroy(&y);CHKERRQ(ierr); 16397929ea7SJunchao Zhang ierr = PetscFinalize(); 16497929ea7SJunchao Zhang return ierr; 16597929ea7SJunchao Zhang } 16697929ea7SJunchao Zhang 16797929ea7SJunchao Zhang /*TEST 16897929ea7SJunchao Zhang 16997929ea7SJunchao Zhang test: 17097929ea7SJunchao Zhang nsize: 4 17197929ea7SJunchao Zhang args: 172*dfd57a17SPierre Jolivet requires: double defined(PETSC_USE_LOG) 17397929ea7SJunchao Zhang 17497929ea7SJunchao Zhang TEST*/ 175