xref: /petsc/src/vec/is/sf/tests/ex15.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
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 
main(int argc,char ** argv)11d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
12d71ae5a4SJacob Faibussowitsch {
1397929ea7SJunchao Zhang   PetscInt        i, n, *ix, *iy, *tomap, start;
1497929ea7SJunchao Zhang   Vec             x, y;
1597929ea7SJunchao Zhang   PetscMPIInt     nproc, rank;
1697929ea7SJunchao Zhang   IS              isx, isy;
1797929ea7SJunchao Zhang   const PetscInt *ranges;
1897929ea7SJunchao Zhang   VecScatter      vscat;
1997929ea7SJunchao Zhang 
2097929ea7SJunchao Zhang   PetscFunctionBegin;
21327415f7SBarry Smith   PetscFunctionBeginUser;
22*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &nproc));
249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
2597929ea7SJunchao Zhang 
2608401ef6SPierre Jolivet   PetscCheck(nproc == 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This test must run with exactly two MPI ranks");
2797929ea7SJunchao Zhang 
2897929ea7SJunchao Zhang   /* ====================================================================
2997929ea7SJunchao Zhang      (1) test VecScatterRemap on a parallel to parallel (PtoP) vecscatter
3097929ea7SJunchao Zhang      ====================================================================
3197929ea7SJunchao Zhang    */
3297929ea7SJunchao Zhang 
3397929ea7SJunchao Zhang   n = 64; /* long enough to trigger memcpy optimizations both in local scatter and remote scatter */
3497929ea7SJunchao Zhang 
3597929ea7SJunchao Zhang   /* create two MPI vectors x, y of length n=64, N=128 */
3677433607SBarry Smith   PetscCall(VecCreateFromOptions(PETSC_COMM_WORLD, NULL, 1, n, PETSC_DECIDE, &x));
379566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &y));
3897929ea7SJunchao Zhang 
3997929ea7SJunchao Zhang   /* Initialize x as {0~127} */
409566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRanges(x, &ranges));
419566063dSJacob Faibussowitsch   for (i = ranges[rank]; i < ranges[rank + 1]; i++) PetscCall(VecSetValue(x, i, (PetscScalar)i, INSERT_VALUES));
429566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(x));
439566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(x));
4497929ea7SJunchao Zhang 
4597929ea7SJunchao 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
4697929ea7SJunchao Zhang      it as general and let PETSc detect the pattern and optimize it. indices in isy are set to make the vecscatter
4797929ea7SJunchao Zhang      have both local scatter and remote scatter (i.e., MPI communication)
4897929ea7SJunchao Zhang    */
499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(n, &ix, n, &iy));
5097929ea7SJunchao Zhang   start = ranges[rank];
5197929ea7SJunchao Zhang   for (i = ranges[rank]; i < ranges[rank + 1]; i++) ix[i - start] = i;
529566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n, ix, PETSC_COPY_VALUES, &isx));
5397929ea7SJunchao Zhang 
549371c9d4SSatish Balay   if (rank == 0) {
559371c9d4SSatish Balay     for (i = 0; i < n; i++) iy[i] = i + 32;
569371c9d4SSatish Balay   } else
579371c9d4SSatish Balay     for (i = 0; i < n / 2; i++) {
589371c9d4SSatish Balay       iy[i]         = i + 96;
599371c9d4SSatish Balay       iy[i + n / 2] = i;
609371c9d4SSatish Balay     }
6197929ea7SJunchao Zhang 
629566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n, iy, PETSC_COPY_VALUES, &isy));
6397929ea7SJunchao Zhang 
6469d47153SPierre Jolivet   /* create a vecscatter that shifts x to the tail by quarter periodically and puts the results in y */
659566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x, isx, y, isy, &vscat));
669566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
679566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
6897929ea7SJunchao Zhang 
6997929ea7SJunchao 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} */
709566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Before VecScatterRemap on PtoP, MPI vector y is:\n"));
719566063dSJacob Faibussowitsch   PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
7297929ea7SJunchao Zhang 
7397929ea7SJunchao Zhang   /* now call the weird subroutine VecScatterRemap to slightly change the vecscatter. It changes where we read vector
7497929ea7SJunchao Zhang      x entries to send out, but does not change the communication pattern (i.e., send/recv pairs and msg lengths).
7597929ea7SJunchao Zhang 
766aad120cSJose 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
7797929ea7SJunchao Zhang      does indices[i] = tomap[indices[i]]. Therefore, after the remap, we read from indices {32~63,0~31} of the local x.
7897929ea7SJunchao Zhang      isy is unchanged. So, we will shift x to {Q2,Q1,Q0,Q3}, that is {64~95,32~63,0~31,96~127}
7997929ea7SJunchao Zhang   */
809566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &tomap));
819371c9d4SSatish Balay   for (i = 0; i < n / 2; i++) {
829371c9d4SSatish Balay     tomap[i]         = i + n / 2;
839371c9d4SSatish Balay     tomap[i + n / 2] = i;
849371c9d4SSatish Balay   };
859566063dSJacob Faibussowitsch   PetscCall(VecScatterRemap(vscat, tomap, NULL));
869566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
879566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
8897929ea7SJunchao Zhang 
8997929ea7SJunchao Zhang   /* view y to check the result. y should be {64~95,32~63,0~31,96~127} */
909566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "After VecScatterRemap on PtoP, MPI vector y is:\n"));
919566063dSJacob Faibussowitsch   PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
9297929ea7SJunchao Zhang 
9397929ea7SJunchao Zhang   /* destroy everything before we recreate them in different types */
949566063dSJacob Faibussowitsch   PetscCall(PetscFree2(ix, iy));
959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
979566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isx));
989566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isy));
999566063dSJacob Faibussowitsch   PetscCall(PetscFree(tomap));
1009566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&vscat));
10197929ea7SJunchao Zhang 
10297929ea7SJunchao Zhang   /* ==========================================================================================
10397929ea7SJunchao Zhang      (2) test VecScatterRemap on a sequential general to sequential general (SGToSG) vecscatter
10497929ea7SJunchao Zhang      ==========================================================================================
10597929ea7SJunchao Zhang    */
10697929ea7SJunchao Zhang   n = 64; /* long enough to trigger memcpy optimizations in local scatter */
10797929ea7SJunchao Zhang 
10897929ea7SJunchao Zhang   /* create two seq vectors x, y of length n */
10977433607SBarry Smith   PetscCall(VecCreateFromOptions(PETSC_COMM_SELF, NULL, 1, n, n, &x));
1109566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &y));
11197929ea7SJunchao Zhang 
11297929ea7SJunchao Zhang   /* Initialize x as {0~63} */
1139566063dSJacob Faibussowitsch   for (i = 0; i < n; i++) PetscCall(VecSetValue(x, i, (PetscScalar)i, INSERT_VALUES));
1149566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(x));
1159566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(x));
11697929ea7SJunchao Zhang 
11797929ea7SJunchao Zhang   /* create two general index sets isx = isy = {0~63}, which are sequential, but we use them as
11897929ea7SJunchao Zhang      general and let PETSc detect the pattern and optimize it */
1199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(n, &ix, n, &iy));
12097929ea7SJunchao Zhang   for (i = 0; i < n; i++) ix[i] = i;
1219566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, ix, PETSC_COPY_VALUES, &isx));
1229566063dSJacob Faibussowitsch   PetscCall(ISDuplicate(isx, &isy));
12397929ea7SJunchao Zhang 
12497929ea7SJunchao Zhang   /* create a vecscatter that just copies x to y */
1259566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x, isx, y, isy, &vscat));
1269566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
1279566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
12897929ea7SJunchao Zhang 
12997929ea7SJunchao Zhang   /* view y to check the result. y should be {0~63} */
1309566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nBefore VecScatterRemap on SGToSG, SEQ vector y is:\n"));
1319566063dSJacob Faibussowitsch   PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
13297929ea7SJunchao Zhang 
13397929ea7SJunchao Zhang   /* now call the weird subroutine VecScatterRemap to slightly change the vecscatter.
13497929ea7SJunchao Zhang 
1356aad120cSJose 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
13697929ea7SJunchao Zhang      does indices[i] = tomap[indices[i]]. Therefore, after the remap, we read from indices{32~63,0~31} of seq x.
13797929ea7SJunchao Zhang    */
1389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &tomap));
1399371c9d4SSatish Balay   for (i = 0; i < n / 2; i++) {
1409371c9d4SSatish Balay     tomap[i]         = i + n / 2;
1419371c9d4SSatish Balay     tomap[i + n / 2] = i;
1429371c9d4SSatish Balay   };
1439566063dSJacob Faibussowitsch   PetscCall(VecScatterRemap(vscat, tomap, NULL));
1449566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
1459566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
14697929ea7SJunchao Zhang 
14797929ea7SJunchao Zhang   /* view y to check the result. y should be {32~63,0~31} */
1489566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "After VecScatterRemap on SGToSG, SEQ vector y is:\n"));
1499566063dSJacob Faibussowitsch   PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
15097929ea7SJunchao Zhang 
15197929ea7SJunchao Zhang   /* destroy everything before we recreate them in different types */
1529566063dSJacob Faibussowitsch   PetscCall(PetscFree2(ix, iy));
1539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
1559566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isx));
1569566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isy));
1579566063dSJacob Faibussowitsch   PetscCall(PetscFree(tomap));
1589566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&vscat));
15997929ea7SJunchao Zhang 
16097929ea7SJunchao Zhang   /* ===================================================================================================
16197929ea7SJunchao Zhang      (3) test VecScatterRemap on a sequential general to sequential stride 1 (SGToSS_Stride1) vecscatter
16297929ea7SJunchao Zhang      ===================================================================================================
16397929ea7SJunchao Zhang    */
16497929ea7SJunchao Zhang   n = 64; /* long enough to trigger memcpy optimizations in local scatter */
16597929ea7SJunchao Zhang 
16697929ea7SJunchao Zhang   /* create two seq vectors x of length n, and y of length n/2 */
16777433607SBarry Smith   PetscCall(VecCreateFromOptions(PETSC_COMM_SELF, NULL, 1, n, n, &x));
16877433607SBarry Smith   PetscCall(VecCreateFromOptions(PETSC_COMM_SELF, NULL, 1, n / 2, n / 2, &y));
16997929ea7SJunchao Zhang 
17097929ea7SJunchao Zhang   /* Initialize x as {0~63} */
1719566063dSJacob Faibussowitsch   for (i = 0; i < n; i++) PetscCall(VecSetValue(x, i, (PetscScalar)i, INSERT_VALUES));
1729566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(x));
1739566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(x));
17497929ea7SJunchao Zhang 
17597929ea7SJunchao Zhang   /* create a general index set isx = {0:63:2}, which actually is a stride IS with first=0, n=32, step=2,
17697929ea7SJunchao Zhang      but we use it as general and let PETSc detect the pattern and optimize it. */
1779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(n / 2, &ix, n / 2, &iy));
17897929ea7SJunchao Zhang   for (i = 0; i < n / 2; i++) ix[i] = i * 2;
1799566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n / 2, ix, PETSC_COPY_VALUES, &isx));
18097929ea7SJunchao Zhang 
18197929ea7SJunchao Zhang   /* create a stride1 index set isy = {0~31}. We intentionally set the step to 1 to trigger optimizations */
1829566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, 32, 0, 1, &isy));
18397929ea7SJunchao Zhang 
18497929ea7SJunchao Zhang   /* create a vecscatter that just copies even entries of x to y */
1859566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x, isx, y, isy, &vscat));
1869566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
1879566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
18897929ea7SJunchao Zhang 
18997929ea7SJunchao Zhang   /* view y to check the result. y should be {0:63:2} */
1909566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nBefore VecScatterRemap on SGToSS_Stride1, SEQ vector y is:\n"));
1919566063dSJacob Faibussowitsch   PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
19297929ea7SJunchao Zhang 
19397929ea7SJunchao Zhang   /* now call the weird subroutine VecScatterRemap to slightly change the vecscatter.
19497929ea7SJunchao Zhang 
1956aad120cSJose 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
19697929ea7SJunchao Zhang      does indices[i] = tomap[indices[i]]. Therefore, after the remap, we read from indices{32:63:2,0:31:2} of seq x.
19797929ea7SJunchao Zhang    */
1989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &tomap));
1999371c9d4SSatish Balay   for (i = 0; i < n / 2; i++) {
2009371c9d4SSatish Balay     tomap[i]         = i + n / 2;
2019371c9d4SSatish Balay     tomap[i + n / 2] = i;
2029371c9d4SSatish Balay   };
2039566063dSJacob Faibussowitsch   PetscCall(VecScatterRemap(vscat, tomap, NULL));
2049566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
2059566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
20697929ea7SJunchao Zhang 
20797929ea7SJunchao Zhang   /* view y to check the result. y should be {32:63:2,0:31:2} */
2089566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "After VecScatterRemap on SGToSS_Stride1, SEQ vector y is:\n"));
2099566063dSJacob Faibussowitsch   PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
21097929ea7SJunchao Zhang 
21197929ea7SJunchao Zhang   /* destroy everything before PetscFinalize */
2129566063dSJacob Faibussowitsch   PetscCall(PetscFree2(ix, iy));
2139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
2159566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isx));
2169566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isy));
2179566063dSJacob Faibussowitsch   PetscCall(PetscFree(tomap));
2189566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&vscat));
21997929ea7SJunchao Zhang 
2209566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
221b122ec5aSJacob Faibussowitsch   return 0;
22297929ea7SJunchao Zhang }
22397929ea7SJunchao Zhang 
22497929ea7SJunchao Zhang /*TEST
22597929ea7SJunchao Zhang 
22697929ea7SJunchao Zhang    test:
22797929ea7SJunchao Zhang       suffix: 1
22897929ea7SJunchao Zhang       nsize: 2
22997929ea7SJunchao Zhang       diff_args: -j
23097929ea7SJunchao Zhang       requires: double
23197929ea7SJunchao Zhang TEST*/
232