xref: /petsc/src/vec/is/sf/tests/ex15.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
197929ea7SJunchao Zhang static char help[] = "  Test VecScatterRemap() on various vecscatter. \n\
297929ea7SJunchao Zhang We may do optimization based on index patterns. After index remapping by VecScatterRemap(), we need to \n\
397929ea7SJunchao Zhang make sure the vecscatter works as expected with the optimizaiton. \n\
497929ea7SJunchao Zhang   VecScatterRemap() does not support all kinds of vecscatters. In addition, it only supports remapping \n\
56aad120cSJose E. Roman entries where we read the data (i.e., todata in parallel scatter, fromdata in sequential scatter). This test \n\
66aad120cSJose E. Roman tests VecScatterRemap on parallel to parallel (PtoP) vecscatter, sequential general to sequential \n\
797929ea7SJunchao Zhang general (SGToSG) vecscatter and sequential general to sequential stride 1 (SGToSS_Stride1) vecscatter.\n\n";
897929ea7SJunchao Zhang 
997929ea7SJunchao Zhang #include <petscvec.h>
1097929ea7SJunchao Zhang 
11*9371c9d4SSatish Balay int main(int argc, char **argv) {
1297929ea7SJunchao Zhang   PetscInt        i, n, *ix, *iy, *tomap, start;
1397929ea7SJunchao Zhang   Vec             x, y;
1497929ea7SJunchao Zhang   PetscMPIInt     nproc, rank;
1597929ea7SJunchao Zhang   IS              isx, isy;
1697929ea7SJunchao Zhang   const PetscInt *ranges;
1797929ea7SJunchao Zhang   VecScatter      vscat;
1897929ea7SJunchao Zhang 
1997929ea7SJunchao Zhang   PetscFunctionBegin;
20327415f7SBarry Smith   PetscFunctionBeginUser;
219566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &nproc));
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
2497929ea7SJunchao Zhang 
2508401ef6SPierre Jolivet   PetscCheck(nproc == 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This test must run with exactly two MPI ranks");
2697929ea7SJunchao Zhang 
2797929ea7SJunchao Zhang   /* ====================================================================
2897929ea7SJunchao Zhang      (1) test VecScatterRemap on a parallel to parallel (PtoP) vecscatter
2997929ea7SJunchao Zhang      ====================================================================
3097929ea7SJunchao Zhang    */
3197929ea7SJunchao Zhang 
3297929ea7SJunchao Zhang   n = 64; /* long enough to trigger memcpy optimizations both in local scatter and remote scatter */
3397929ea7SJunchao Zhang 
3497929ea7SJunchao Zhang   /* create two MPI vectors x, y of length n=64, N=128 */
359566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PETSC_COMM_WORLD, n, PETSC_DECIDE, &x));
369566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &y));
3797929ea7SJunchao Zhang 
3897929ea7SJunchao Zhang   /* Initialize x as {0~127} */
399566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRanges(x, &ranges));
409566063dSJacob Faibussowitsch   for (i = ranges[rank]; i < ranges[rank + 1]; i++) PetscCall(VecSetValue(x, i, (PetscScalar)i, INSERT_VALUES));
419566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(x));
429566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(x));
4397929ea7SJunchao Zhang 
4497929ea7SJunchao Zhang   /* create two general index sets isx = {0~127} and isy = {32~63,64~95,96~127,0~31}. isx is sequential, but we use
4597929ea7SJunchao Zhang      it as general and let PETSc detect the pattern and optimize it. indices in isy are set to make the vecscatter
4697929ea7SJunchao Zhang      have both local scatter and remote scatter (i.e., MPI communication)
4797929ea7SJunchao Zhang    */
489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(n, &ix, n, &iy));
4997929ea7SJunchao Zhang   start = ranges[rank];
5097929ea7SJunchao Zhang   for (i = ranges[rank]; i < ranges[rank + 1]; i++) ix[i - start] = i;
519566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n, ix, PETSC_COPY_VALUES, &isx));
5297929ea7SJunchao Zhang 
53*9371c9d4SSatish Balay   if (rank == 0) {
54*9371c9d4SSatish Balay     for (i = 0; i < n; i++) iy[i] = i + 32;
55*9371c9d4SSatish Balay   } else
56*9371c9d4SSatish Balay     for (i = 0; i < n / 2; i++) {
57*9371c9d4SSatish Balay       iy[i]         = i + 96;
58*9371c9d4SSatish Balay       iy[i + n / 2] = i;
59*9371c9d4SSatish Balay     }
6097929ea7SJunchao Zhang 
619566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n, iy, PETSC_COPY_VALUES, &isy));
6297929ea7SJunchao Zhang 
6397929ea7SJunchao Zhang   /* create a vecscatter that shifts x to the tail by quater periodically and puts the results in y */
649566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x, isx, y, isy, &vscat));
659566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
669566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
6797929ea7SJunchao Zhang 
6897929ea7SJunchao Zhang   /* view y to check the result. y should be {Q3,Q0,Q1,Q2} of x, that is {96~127,0~31,32~63,64~95} */
699566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Before VecScatterRemap on PtoP, MPI vector y is:\n"));
709566063dSJacob Faibussowitsch   PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
7197929ea7SJunchao Zhang 
7297929ea7SJunchao Zhang   /* now call the weird subroutine VecScatterRemap to slightly change the vecscatter. It changes where we read vector
7397929ea7SJunchao Zhang      x entries to send out, but does not change the communication pattern (i.e., send/recv pairs and msg lengths).
7497929ea7SJunchao Zhang 
756aad120cSJose E. Roman      We create tomap as {32~63,0~31}. Originally, we read from indices {0~64} of the local x to send out. The remap
7697929ea7SJunchao Zhang      does indices[i] = tomap[indices[i]]. Therefore, after the remap, we read from indices {32~63,0~31} of the local x.
7797929ea7SJunchao Zhang      isy is unchanged. So, we will shift x to {Q2,Q1,Q0,Q3}, that is {64~95,32~63,0~31,96~127}
7897929ea7SJunchao Zhang   */
799566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &tomap));
80*9371c9d4SSatish Balay   for (i = 0; i < n / 2; i++) {
81*9371c9d4SSatish Balay     tomap[i]         = i + n / 2;
82*9371c9d4SSatish Balay     tomap[i + n / 2] = i;
83*9371c9d4SSatish Balay   };
849566063dSJacob Faibussowitsch   PetscCall(VecScatterRemap(vscat, tomap, NULL));
859566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
869566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
8797929ea7SJunchao Zhang 
8897929ea7SJunchao Zhang   /* view y to check the result. y should be {64~95,32~63,0~31,96~127} */
899566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "After VecScatterRemap on PtoP, MPI vector y is:\n"));
909566063dSJacob Faibussowitsch   PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
9197929ea7SJunchao Zhang 
9297929ea7SJunchao Zhang   /* destroy everything before we recreate them in different types */
939566063dSJacob Faibussowitsch   PetscCall(PetscFree2(ix, iy));
949566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
969566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isx));
979566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isy));
989566063dSJacob Faibussowitsch   PetscCall(PetscFree(tomap));
999566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&vscat));
10097929ea7SJunchao Zhang 
10197929ea7SJunchao Zhang   /* ==========================================================================================
10297929ea7SJunchao Zhang      (2) test VecScatterRemap on a sequential general to sequential general (SGToSG) vecscatter
10397929ea7SJunchao Zhang      ==========================================================================================
10497929ea7SJunchao Zhang    */
10597929ea7SJunchao Zhang   n = 64; /* long enough to trigger memcpy optimizations in local scatter */
10697929ea7SJunchao Zhang 
10797929ea7SJunchao Zhang   /* create two seq vectors x, y of length n */
1089566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &x));
1099566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &y));
11097929ea7SJunchao Zhang 
11197929ea7SJunchao Zhang   /* Initialize x as {0~63} */
1129566063dSJacob Faibussowitsch   for (i = 0; i < n; i++) PetscCall(VecSetValue(x, i, (PetscScalar)i, INSERT_VALUES));
1139566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(x));
1149566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(x));
11597929ea7SJunchao Zhang 
11697929ea7SJunchao Zhang   /* create two general index sets isx = isy = {0~63}, which are sequential, but we use them as
11797929ea7SJunchao Zhang      general and let PETSc detect the pattern and optimize it */
1189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(n, &ix, n, &iy));
11997929ea7SJunchao Zhang   for (i = 0; i < n; i++) ix[i] = i;
1209566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, ix, PETSC_COPY_VALUES, &isx));
1219566063dSJacob Faibussowitsch   PetscCall(ISDuplicate(isx, &isy));
12297929ea7SJunchao Zhang 
12397929ea7SJunchao Zhang   /* create a vecscatter that just copies x to y */
1249566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x, isx, y, isy, &vscat));
1259566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
1269566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
12797929ea7SJunchao Zhang 
12897929ea7SJunchao Zhang   /* view y to check the result. y should be {0~63} */
1299566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nBefore VecScatterRemap on SGToSG, SEQ vector y is:\n"));
1309566063dSJacob Faibussowitsch   PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
13197929ea7SJunchao Zhang 
13297929ea7SJunchao Zhang   /* now call the weird subroutine VecScatterRemap to slightly change the vecscatter.
13397929ea7SJunchao Zhang 
1346aad120cSJose E. Roman      Create tomap as {32~63,0~31}. Originally, we read from indices {0~64} of seq x to write to y. The remap
13597929ea7SJunchao Zhang      does indices[i] = tomap[indices[i]]. Therefore, after the remap, we read from indices{32~63,0~31} of seq x.
13697929ea7SJunchao Zhang    */
1379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &tomap));
138*9371c9d4SSatish Balay   for (i = 0; i < n / 2; i++) {
139*9371c9d4SSatish Balay     tomap[i]         = i + n / 2;
140*9371c9d4SSatish Balay     tomap[i + n / 2] = i;
141*9371c9d4SSatish Balay   };
1429566063dSJacob Faibussowitsch   PetscCall(VecScatterRemap(vscat, tomap, NULL));
1439566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
1449566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
14597929ea7SJunchao Zhang 
14697929ea7SJunchao Zhang   /* view y to check the result. y should be {32~63,0~31} */
1479566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "After VecScatterRemap on SGToSG, SEQ vector y is:\n"));
1489566063dSJacob Faibussowitsch   PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
14997929ea7SJunchao Zhang 
15097929ea7SJunchao Zhang   /* destroy everything before we recreate them in different types */
1519566063dSJacob Faibussowitsch   PetscCall(PetscFree2(ix, iy));
1529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
1549566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isx));
1559566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isy));
1569566063dSJacob Faibussowitsch   PetscCall(PetscFree(tomap));
1579566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&vscat));
15897929ea7SJunchao Zhang 
15997929ea7SJunchao Zhang   /* ===================================================================================================
16097929ea7SJunchao Zhang      (3) test VecScatterRemap on a sequential general to sequential stride 1 (SGToSS_Stride1) vecscatter
16197929ea7SJunchao Zhang      ===================================================================================================
16297929ea7SJunchao Zhang    */
16397929ea7SJunchao Zhang   n = 64; /* long enough to trigger memcpy optimizations in local scatter */
16497929ea7SJunchao Zhang 
16597929ea7SJunchao Zhang   /* create two seq vectors x of length n, and y of length n/2 */
1669566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &x));
1679566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n / 2, &y));
16897929ea7SJunchao Zhang 
16997929ea7SJunchao Zhang   /* Initialize x as {0~63} */
1709566063dSJacob Faibussowitsch   for (i = 0; i < n; i++) PetscCall(VecSetValue(x, i, (PetscScalar)i, INSERT_VALUES));
1719566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(x));
1729566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(x));
17397929ea7SJunchao Zhang 
17497929ea7SJunchao Zhang   /* create a general index set isx = {0:63:2}, which actually is a stride IS with first=0, n=32, step=2,
17597929ea7SJunchao Zhang      but we use it as general and let PETSc detect the pattern and optimize it. */
1769566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(n / 2, &ix, n / 2, &iy));
17797929ea7SJunchao Zhang   for (i = 0; i < n / 2; i++) ix[i] = i * 2;
1789566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n / 2, ix, PETSC_COPY_VALUES, &isx));
17997929ea7SJunchao Zhang 
18097929ea7SJunchao Zhang   /* create a stride1 index set isy = {0~31}. We intentionally set the step to 1 to trigger optimizations */
1819566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, 32, 0, 1, &isy));
18297929ea7SJunchao Zhang 
18397929ea7SJunchao Zhang   /* create a vecscatter that just copies even entries of x to y */
1849566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x, isx, y, isy, &vscat));
1859566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
1869566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
18797929ea7SJunchao Zhang 
18897929ea7SJunchao Zhang   /* view y to check the result. y should be {0:63:2} */
1899566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nBefore VecScatterRemap on SGToSS_Stride1, SEQ vector y is:\n"));
1909566063dSJacob Faibussowitsch   PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
19197929ea7SJunchao Zhang 
19297929ea7SJunchao Zhang   /* now call the weird subroutine VecScatterRemap to slightly change the vecscatter.
19397929ea7SJunchao Zhang 
1946aad120cSJose E. Roman      Create tomap as {32~63,0~31}. Originally, we read from indices{0:63:2} of seq x to write to y. The remap
19597929ea7SJunchao Zhang      does indices[i] = tomap[indices[i]]. Therefore, after the remap, we read from indices{32:63:2,0:31:2} of seq x.
19697929ea7SJunchao Zhang    */
1979566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &tomap));
198*9371c9d4SSatish Balay   for (i = 0; i < n / 2; i++) {
199*9371c9d4SSatish Balay     tomap[i]         = i + n / 2;
200*9371c9d4SSatish Balay     tomap[i + n / 2] = i;
201*9371c9d4SSatish Balay   };
2029566063dSJacob Faibussowitsch   PetscCall(VecScatterRemap(vscat, tomap, NULL));
2039566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
2049566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
20597929ea7SJunchao Zhang 
20697929ea7SJunchao Zhang   /* view y to check the result. y should be {32:63:2,0:31:2} */
2079566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "After VecScatterRemap on SGToSS_Stride1, SEQ vector y is:\n"));
2089566063dSJacob Faibussowitsch   PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
20997929ea7SJunchao Zhang 
21097929ea7SJunchao Zhang   /* destroy everything before PetscFinalize */
2119566063dSJacob Faibussowitsch   PetscCall(PetscFree2(ix, iy));
2129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
2149566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isx));
2159566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isy));
2169566063dSJacob Faibussowitsch   PetscCall(PetscFree(tomap));
2179566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&vscat));
21897929ea7SJunchao Zhang 
2199566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
220b122ec5aSJacob Faibussowitsch   return 0;
22197929ea7SJunchao Zhang }
22297929ea7SJunchao Zhang 
22397929ea7SJunchao Zhang /*TEST
22497929ea7SJunchao Zhang 
22597929ea7SJunchao Zhang    test:
22697929ea7SJunchao Zhang       suffix: 1
22797929ea7SJunchao Zhang       nsize: 2
22897929ea7SJunchao Zhang       diff_args: -j
22997929ea7SJunchao Zhang       requires: double
23097929ea7SJunchao Zhang TEST*/
231