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