1*97929ea7SJunchao Zhang 2*97929ea7SJunchao Zhang static char help[]= "Test event log of VecScatter with various block sizes\n\n"; 3*97929ea7SJunchao Zhang 4*97929ea7SJunchao Zhang #include <petscvec.h> 5*97929ea7SJunchao Zhang 6*97929ea7SJunchao Zhang int main(int argc,char **argv) 7*97929ea7SJunchao Zhang { 8*97929ea7SJunchao Zhang PetscErrorCode ierr; 9*97929ea7SJunchao Zhang PetscInt i,j,low,high,n=256,N,errors,tot_errors; 10*97929ea7SJunchao Zhang PetscInt bs=1,ix[2],iy[2]; 11*97929ea7SJunchao Zhang PetscMPIInt nproc,rank; 12*97929ea7SJunchao Zhang PetscScalar *xval; 13*97929ea7SJunchao Zhang const PetscScalar *yval; 14*97929ea7SJunchao Zhang Vec x,y; 15*97929ea7SJunchao Zhang IS isx,isy; 16*97929ea7SJunchao Zhang VecScatter ctx; 17*97929ea7SJunchao Zhang const PetscInt niter = 10; 18*97929ea7SJunchao Zhang #if defined(PETSC_USE_LOG) 19*97929ea7SJunchao Zhang PetscLogStage stage1,stage2; 20*97929ea7SJunchao Zhang PetscLogEvent event1,event2; 21*97929ea7SJunchao Zhang PetscLogDouble numMessages,messageLength; 22*97929ea7SJunchao Zhang PetscEventPerfInfo eventInfo; 23*97929ea7SJunchao Zhang PetscInt tot_msg,tot_len,avg_len; 24*97929ea7SJunchao Zhang #endif 25*97929ea7SJunchao Zhang 26*97929ea7SJunchao Zhang PetscFunctionBegin; 27*97929ea7SJunchao Zhang ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 28*97929ea7SJunchao Zhang ierr = PetscLogDefaultBegin();CHKERRQ(ierr); 29*97929ea7SJunchao Zhang ierr = MPI_Comm_size(PETSC_COMM_WORLD,&nproc);CHKERRQ(ierr); 30*97929ea7SJunchao Zhang ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 31*97929ea7SJunchao Zhang 32*97929ea7SJunchao Zhang ierr = PetscLogStageRegister("Scatter(bs=1)", &stage1);CHKERRQ(ierr); 33*97929ea7SJunchao Zhang ierr = PetscLogEventRegister("VecScatter(bs=1)", PETSC_OBJECT_CLASSID, &event1);CHKERRQ(ierr); 34*97929ea7SJunchao Zhang ierr = PetscLogStageRegister("Scatter(bs=4)", &stage2);CHKERRQ(ierr); 35*97929ea7SJunchao Zhang ierr = PetscLogEventRegister("VecScatter(bs=4)", PETSC_OBJECT_CLASSID, &event2);CHKERRQ(ierr); 36*97929ea7SJunchao Zhang 37*97929ea7SJunchao Zhang /* Create a parallel vector x and a sequential vector y */ 38*97929ea7SJunchao Zhang ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 39*97929ea7SJunchao Zhang ierr = VecSetSizes(x,n,PETSC_DECIDE);CHKERRQ(ierr); 40*97929ea7SJunchao Zhang ierr = VecSetFromOptions(x);CHKERRQ(ierr); 41*97929ea7SJunchao Zhang ierr = VecGetOwnershipRange(x,&low,&high);CHKERRQ(ierr); 42*97929ea7SJunchao Zhang ierr = VecGetSize(x,&N);CHKERRQ(ierr); 43*97929ea7SJunchao Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,n,&y);CHKERRQ(ierr); 44*97929ea7SJunchao Zhang 45*97929ea7SJunchao Zhang /*======================================= 46*97929ea7SJunchao Zhang test VecScatter with bs = 1 47*97929ea7SJunchao Zhang ======================================*/ 48*97929ea7SJunchao Zhang 49*97929ea7SJunchao Zhang /* the code works as if we are going to do 3-point stencil computations on a 1D domain x, 50*97929ea7SJunchao Zhang which has periodic boundary conditions but the two ghosts are scatterred to beginning of y. 51*97929ea7SJunchao Zhang */ 52*97929ea7SJunchao Zhang bs = 1; 53*97929ea7SJunchao Zhang ix[0] = rank? low-1 : N-1; /* ix[] contains global indices of the two ghost points */ 54*97929ea7SJunchao Zhang ix[1] = (rank != nproc-1)? high : 0; 55*97929ea7SJunchao Zhang iy[0] = 0; 56*97929ea7SJunchao Zhang iy[1] = 1; 57*97929ea7SJunchao Zhang 58*97929ea7SJunchao Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,2,ix,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr); 59*97929ea7SJunchao Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,2,iy,PETSC_COPY_VALUES,&isy);CHKERRQ(ierr); 60*97929ea7SJunchao Zhang ierr = VecScatterCreate(x,isx,y,isy,&ctx);CHKERRQ(ierr); 61*97929ea7SJunchao Zhang ierr = VecScatterSetUp(ctx);CHKERRQ(ierr); 62*97929ea7SJunchao Zhang 63*97929ea7SJunchao Zhang ierr = PetscLogStagePush(stage1);CHKERRQ(ierr); 64*97929ea7SJunchao Zhang ierr = PetscLogEventBegin(event1,0,0,0,0);CHKERRQ(ierr); 65*97929ea7SJunchao Zhang errors = 0; 66*97929ea7SJunchao Zhang for (i=0; i<niter; i++) { 67*97929ea7SJunchao Zhang /* set x = 0+i, 1+i, 2+i, ..., N-1+i */ 68*97929ea7SJunchao Zhang ierr = VecGetArray(x,&xval);CHKERRQ(ierr); 69*97929ea7SJunchao Zhang for (j=0; j<n; j++) xval[j] = (PetscScalar)(low+j+i); 70*97929ea7SJunchao Zhang ierr = VecRestoreArray(x,&xval);CHKERRQ(ierr); 71*97929ea7SJunchao Zhang /* scatter the ghosts to y */ 72*97929ea7SJunchao Zhang ierr = VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 73*97929ea7SJunchao Zhang ierr = VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 74*97929ea7SJunchao Zhang /* check if y has correct values */ 75*97929ea7SJunchao Zhang ierr = VecGetArrayRead(y,&yval);CHKERRQ(ierr); 76*97929ea7SJunchao Zhang if ((PetscInt)PetscRealPart(yval[0]) != ix[0]+i) errors++; 77*97929ea7SJunchao Zhang if ((PetscInt)PetscRealPart(yval[1]) != ix[1]+i) errors++; 78*97929ea7SJunchao Zhang ierr = VecRestoreArrayRead(y,&yval);CHKERRQ(ierr); 79*97929ea7SJunchao Zhang } 80*97929ea7SJunchao Zhang ierr = PetscLogEventEnd(event1,0,0,0,0);CHKERRQ(ierr); 81*97929ea7SJunchao Zhang ierr = PetscLogStagePop();CHKERRQ(ierr); 82*97929ea7SJunchao Zhang 83*97929ea7SJunchao Zhang /* check if we found wrong values on any processors */ 84*97929ea7SJunchao Zhang ierr = MPI_Allreduce(&errors,&tot_errors,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);CHKERRQ(ierr); 85*97929ea7SJunchao Zhang if (tot_errors) { ierr = PetscPrintf(PETSC_COMM_WORLD,"Error: wrong values were scatterred in vecscatter with bs = %d\n",bs);CHKERRQ(ierr); } 86*97929ea7SJunchao Zhang 87*97929ea7SJunchao Zhang /* print out event log of VecScatter(bs=1) */ 88*97929ea7SJunchao Zhang #if defined(PETSC_USE_LOG) 89*97929ea7SJunchao Zhang ierr = PetscLogEventGetPerfInfo(stage1,event1,&eventInfo);CHKERRQ(ierr); 90*97929ea7SJunchao Zhang ierr = MPI_Allreduce(&eventInfo.numMessages, &numMessages, 1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD);CHKERRQ(ierr); 91*97929ea7SJunchao Zhang ierr = MPI_Allreduce(&eventInfo.messageLength,&messageLength,1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD);CHKERRQ(ierr); 92*97929ea7SJunchao Zhang tot_msg = (PetscInt)numMessages*0.5; /* two MPI calls (Send & Recv) per message */ 93*97929ea7SJunchao Zhang tot_len = (PetscInt)messageLength*0.5; 94*97929ea7SJunchao Zhang avg_len = tot_msg? (PetscInt)(messageLength/numMessages) : 0; 95*97929ea7SJunchao Zhang /* when nproc > 2, tot_msg = 2*nproc*niter, tot_len = tot_msg*sizeof(PetscScalar)*bs */ 96*97929ea7SJunchao 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); 97*97929ea7SJunchao Zhang #endif 98*97929ea7SJunchao Zhang 99*97929ea7SJunchao Zhang ierr = ISDestroy(&isx);CHKERRQ(ierr); 100*97929ea7SJunchao Zhang ierr = ISDestroy(&isy);CHKERRQ(ierr); 101*97929ea7SJunchao Zhang ierr = VecScatterDestroy(&ctx);CHKERRQ(ierr); 102*97929ea7SJunchao Zhang 103*97929ea7SJunchao Zhang /*======================================= 104*97929ea7SJunchao Zhang test VecScatter with bs = 4 105*97929ea7SJunchao Zhang ======================================*/ 106*97929ea7SJunchao Zhang 107*97929ea7SJunchao Zhang /* similar to the 3-point stencil above, except that this time a ghost is a block */ 108*97929ea7SJunchao Zhang bs = 4; /* n needs to be a multiple of bs to make the following code work */ 109*97929ea7SJunchao Zhang ix[0] = rank? low/bs-1 : N/bs-1; /* ix[] contains global indices of the two ghost blocks */ 110*97929ea7SJunchao Zhang ix[1] = (rank != nproc-1)? high/bs : 0; 111*97929ea7SJunchao Zhang iy[0] = 0; 112*97929ea7SJunchao Zhang iy[1] = 1; 113*97929ea7SJunchao Zhang 114*97929ea7SJunchao Zhang ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2,ix,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr); 115*97929ea7SJunchao Zhang ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2,iy,PETSC_COPY_VALUES,&isy);CHKERRQ(ierr); 116*97929ea7SJunchao Zhang 117*97929ea7SJunchao Zhang ierr = VecScatterCreate(x,isx,y,isy,&ctx);CHKERRQ(ierr); 118*97929ea7SJunchao Zhang /* Call SetUp explicitly, otherwise messages in implicit SetUp will be counted in events below */ 119*97929ea7SJunchao Zhang ierr = VecScatterSetUp(ctx);CHKERRQ(ierr); 120*97929ea7SJunchao Zhang 121*97929ea7SJunchao Zhang ierr = PetscLogStagePush(stage2);CHKERRQ(ierr); 122*97929ea7SJunchao Zhang ierr = PetscLogEventBegin(event2,0,0,0,0);CHKERRQ(ierr); 123*97929ea7SJunchao Zhang errors = 0; 124*97929ea7SJunchao Zhang for (i=0; i<niter; i++) { 125*97929ea7SJunchao Zhang /* set x = 0+i, 1+i, 2+i, ..., N-1+i */ 126*97929ea7SJunchao Zhang ierr = VecGetArray(x,&xval);CHKERRQ(ierr); 127*97929ea7SJunchao Zhang for (j=0; j<n; j++) xval[j] = (PetscScalar)(low+j+i); 128*97929ea7SJunchao Zhang ierr = VecRestoreArray(x,&xval);CHKERRQ(ierr); 129*97929ea7SJunchao Zhang /* scatter the ghost blocks to y */ 130*97929ea7SJunchao Zhang ierr = VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 131*97929ea7SJunchao Zhang ierr = VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 132*97929ea7SJunchao Zhang /* check if y has correct values */ 133*97929ea7SJunchao Zhang ierr = VecGetArrayRead(y,&yval);CHKERRQ(ierr); 134*97929ea7SJunchao Zhang if ((PetscInt)PetscRealPart(yval[0]) != ix[0]*bs+i) errors++; 135*97929ea7SJunchao Zhang if ((PetscInt)PetscRealPart(yval[bs]) != ix[1]*bs+i) errors++; 136*97929ea7SJunchao Zhang ierr = VecRestoreArrayRead(y,&yval);CHKERRQ(ierr); 137*97929ea7SJunchao Zhang } 138*97929ea7SJunchao Zhang ierr = PetscLogEventEnd(event2,0,0,0,0);CHKERRQ(ierr); 139*97929ea7SJunchao Zhang ierr = PetscLogStagePop();CHKERRQ(ierr); 140*97929ea7SJunchao Zhang 141*97929ea7SJunchao Zhang /* check if we found wrong values on any processors */ 142*97929ea7SJunchao Zhang ierr = MPI_Allreduce(&errors,&tot_errors,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);CHKERRQ(ierr); 143*97929ea7SJunchao Zhang if (tot_errors) { ierr = PetscPrintf(PETSC_COMM_WORLD,"Error: wrong values were scatterred in vecscatter with bs = %d\n",bs);CHKERRQ(ierr); } 144*97929ea7SJunchao Zhang 145*97929ea7SJunchao Zhang /* print out event log of VecScatter(bs=4) */ 146*97929ea7SJunchao Zhang #if defined(PETSC_USE_LOG) 147*97929ea7SJunchao Zhang ierr = PetscLogEventGetPerfInfo(stage2,event2,&eventInfo);CHKERRQ(ierr); 148*97929ea7SJunchao Zhang ierr = MPI_Allreduce(&eventInfo.numMessages, &numMessages, 1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD);CHKERRQ(ierr); 149*97929ea7SJunchao Zhang ierr = MPI_Allreduce(&eventInfo.messageLength,&messageLength,1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD);CHKERRQ(ierr); 150*97929ea7SJunchao Zhang tot_msg = (PetscInt)numMessages*0.5; /* two MPI calls (Send & Recv) per message */ 151*97929ea7SJunchao Zhang tot_len = (PetscInt)messageLength*0.5; 152*97929ea7SJunchao Zhang avg_len = tot_msg? (PetscInt)(messageLength/numMessages) : 0; 153*97929ea7SJunchao Zhang /* when nproc > 2, tot_msg = 2*nproc*niter, tot_len = tot_msg*sizeof(PetscScalar)*bs */ 154*97929ea7SJunchao 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); 155*97929ea7SJunchao Zhang #endif 156*97929ea7SJunchao Zhang 157*97929ea7SJunchao Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"Program finished\n");CHKERRQ(ierr); 158*97929ea7SJunchao Zhang ierr = ISDestroy(&isx);CHKERRQ(ierr); 159*97929ea7SJunchao Zhang ierr = ISDestroy(&isy);CHKERRQ(ierr); 160*97929ea7SJunchao Zhang ierr = VecScatterDestroy(&ctx);CHKERRQ(ierr); 161*97929ea7SJunchao Zhang ierr = VecDestroy(&x);CHKERRQ(ierr); 162*97929ea7SJunchao Zhang ierr = VecDestroy(&y);CHKERRQ(ierr); 163*97929ea7SJunchao Zhang ierr = PetscFinalize(); 164*97929ea7SJunchao Zhang return ierr; 165*97929ea7SJunchao Zhang } 166*97929ea7SJunchao Zhang 167*97929ea7SJunchao Zhang /*TEST 168*97929ea7SJunchao Zhang 169*97929ea7SJunchao Zhang test: 170*97929ea7SJunchao Zhang nsize: 4 171*97929ea7SJunchao Zhang args: 172*97929ea7SJunchao Zhang requires: double define(PETSC_USE_LOG) 173*97929ea7SJunchao Zhang 174*97929ea7SJunchao Zhang TEST*/ 175