xref: /petsc/src/vec/is/sf/tests/ex14.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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;
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogDefaultBegin());
29*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&nproc));
30*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
3197929ea7SJunchao Zhang 
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStageRegister("Scatter(bs=1)", &stage1));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventRegister("VecScatter(bs=1)", PETSC_OBJECT_CLASSID, &event1));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStageRegister("Scatter(bs=4)", &stage2));
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventRegister("VecScatter(bs=4)", PETSC_OBJECT_CLASSID, &event2));
3697929ea7SJunchao Zhang 
3797929ea7SJunchao Zhang   /* Create a parallel vector x and a sequential vector y */
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(x,n,PETSC_DECIDE));
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(x));
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(x,&low,&high));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(x,&N));
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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 
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,2,ix,PETSC_COPY_VALUES,&isx));
59*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,2,iy,PETSC_COPY_VALUES,&isy));
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(x,isx,y,isy,&ctx));
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterSetUp(ctx));
6297929ea7SJunchao Zhang 
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePush(stage1));
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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 */
68*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(x,&xval));
6997929ea7SJunchao Zhang     for (j=0; j<n; j++) xval[j] = (PetscScalar)(low+j+i);
70*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(x,&xval));
7197929ea7SJunchao Zhang     /* scatter the ghosts to y */
72*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD));
73*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD));
7497929ea7SJunchao Zhang     /* check if y has correct values */
75*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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++;
78*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(y,&yval));
7997929ea7SJunchao Zhang   }
80*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(event1,0,0,0,0));
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePop());
8297929ea7SJunchao Zhang 
8397929ea7SJunchao Zhang   /* check if we found wrong values on any processors */
84*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allreduce(&errors,&tot_errors,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD));
85*5f80ce2aSJacob Faibussowitsch   if (tot_errors) CHKERRQ(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)
89*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventGetPerfInfo(stage1,event1,&eventInfo));
90*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allreduce(&eventInfo.numMessages,  &numMessages,  1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD));
91*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(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 */
96*5f80ce2aSJacob 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));
9797929ea7SJunchao Zhang #endif
9897929ea7SJunchao Zhang 
99*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isx));
100*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isy));
101*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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 
114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,2,ix,PETSC_COPY_VALUES,&isx));
115*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,2,iy,PETSC_COPY_VALUES,&isy));
11697929ea7SJunchao Zhang 
117*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(x,isx,y,isy,&ctx));
11897929ea7SJunchao Zhang    /* Call SetUp explicitly, otherwise messages in implicit SetUp will be counted in events below */
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterSetUp(ctx));
12097929ea7SJunchao Zhang 
121*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePush(stage2));
122*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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 */
126*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(x,&xval));
12797929ea7SJunchao Zhang     for (j=0; j<n; j++) xval[j] = (PetscScalar)(low+j+i);
128*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(x,&xval));
12997929ea7SJunchao Zhang     /* scatter the ghost blocks to y */
130*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD));
131*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD));
13297929ea7SJunchao Zhang     /* check if y has correct values */
133*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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++;
136*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(y,&yval));
13797929ea7SJunchao Zhang   }
138*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(event2,0,0,0,0));
139*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePop());
14097929ea7SJunchao Zhang 
14197929ea7SJunchao Zhang   /* check if we found wrong values on any processors */
142*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allreduce(&errors,&tot_errors,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD));
143*5f80ce2aSJacob Faibussowitsch   if (tot_errors) CHKERRQ(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)
147*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventGetPerfInfo(stage2,event2,&eventInfo));
148*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allreduce(&eventInfo.numMessages,  &numMessages,  1,MPIU_PETSCLOGDOUBLE,MPI_SUM,PETSC_COMM_WORLD));
149*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(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 */
154*5f80ce2aSJacob 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));
15597929ea7SJunchao Zhang #endif
15697929ea7SJunchao Zhang 
157*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Program finished\n"));
158*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isx));
159*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isy));
160*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&ctx));
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
162*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
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:
172dfd57a17SPierre Jolivet       requires: double defined(PETSC_USE_LOG)
17397929ea7SJunchao Zhang 
17497929ea7SJunchao Zhang TEST*/
175