xref: /petsc/src/vec/is/sf/tests/ex15.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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\
597929ea7SJunchao Zhang entries where we read the data (i.e., todata in paralle scatter, fromdata in sequential scatter). This test \n\
697929ea7SJunchao Zhang tests VecScatterRemap on parallel to paralle (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 
1197929ea7SJunchao Zhang int main(int argc,char **argv)
1297929ea7SJunchao Zhang {
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;
21*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
225f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&nproc));
235f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
2497929ea7SJunchao Zhang 
252c71b3e2SJacob Faibussowitsch   PetscCheckFalse(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 */
355f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,n,PETSC_DECIDE,&x));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&y));
3797929ea7SJunchao Zhang 
3897929ea7SJunchao Zhang   /* Initialize x as {0~127} */
395f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRanges(x,&ranges));
405f80ce2aSJacob Faibussowitsch   for (i=ranges[rank]; i<ranges[rank+1]; i++) CHKERRQ(VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(x));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(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    */
485f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
515f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,n,ix,PETSC_COPY_VALUES,&isx));
5297929ea7SJunchao Zhang 
53dd400576SPatrick Sanan   if (rank == 0) { for (i=0; i<n; i++) iy[i] = i+32; }
5497929ea7SJunchao Zhang   else for (i=0; i<n/2; i++) { iy[i] = i+96; iy[i+n/2] = i; }
5597929ea7SJunchao Zhang 
565f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,n,iy,PETSC_COPY_VALUES,&isy));
5797929ea7SJunchao Zhang 
5897929ea7SJunchao Zhang   /* create a vecscatter that shifts x to the tail by quater periodically and puts the results in y */
595f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(x,isx,y,isy,&vscat));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
6297929ea7SJunchao Zhang 
6397929ea7SJunchao 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} */
645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Before VecScatterRemap on PtoP, MPI vector y is:\n"));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
6697929ea7SJunchao Zhang 
6797929ea7SJunchao Zhang   /* now call the weird subroutine VecScatterRemap to slightly change the vecscatter. It changes where we read vector
6897929ea7SJunchao Zhang      x entries to send out, but does not change the communication pattern (i.e., send/recv pairs and msg lengths).
6997929ea7SJunchao Zhang 
7097929ea7SJunchao Zhang      We create tomap as {32~63,0~31}. Originaly, we read from indices {0~64} of the local x to send out. The remap
7197929ea7SJunchao Zhang      does indices[i] = tomap[indices[i]]. Therefore, after the remap, we read from indices {32~63,0~31} of the local x.
7297929ea7SJunchao Zhang      isy is unchanged. So, we will shift x to {Q2,Q1,Q0,Q3}, that is {64~95,32~63,0~31,96~127}
7397929ea7SJunchao Zhang   */
745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&tomap));
7597929ea7SJunchao Zhang   for (i=0; i<n/2; i++) { tomap[i] = i+n/2; tomap[i+n/2] = i; };
765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterRemap(vscat,tomap,NULL));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
7997929ea7SJunchao Zhang 
8097929ea7SJunchao Zhang   /* view y to check the result. y should be {64~95,32~63,0~31,96~127} */
815f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"After VecScatterRemap on PtoP, MPI vector y is:\n"));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
8397929ea7SJunchao Zhang 
8497929ea7SJunchao Zhang   /* destroy everything before we recreate them in different types */
855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(ix,iy));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isx));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isy));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(tomap));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&vscat));
9297929ea7SJunchao Zhang 
9397929ea7SJunchao Zhang   /* ==========================================================================================
9497929ea7SJunchao Zhang      (2) test VecScatterRemap on a sequential general to sequential general (SGToSG) vecscatter
9597929ea7SJunchao Zhang      ==========================================================================================
9697929ea7SJunchao Zhang    */
9797929ea7SJunchao Zhang   n = 64; /* long enough to trigger memcpy optimizations in local scatter */
9897929ea7SJunchao Zhang 
9997929ea7SJunchao Zhang   /* create two seq vectors x, y of length n */
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&x));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&y));
10297929ea7SJunchao Zhang 
10397929ea7SJunchao Zhang   /* Initialize x as {0~63} */
1045f80ce2aSJacob Faibussowitsch   for (i=0; i<n; i++) CHKERRQ(VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(x));
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(x));
10797929ea7SJunchao Zhang 
10897929ea7SJunchao Zhang   /* create two general index sets isx = isy = {0~63}, which are sequential, but we use them as
10997929ea7SJunchao Zhang      general and let PETSc detect the pattern and optimize it */
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(n,&ix,n,&iy));
11197929ea7SJunchao Zhang   for (i=0; i<n; i++) ix[i] = i;
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,n,ix,PETSC_COPY_VALUES,&isx));
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDuplicate(isx,&isy));
11497929ea7SJunchao Zhang 
11597929ea7SJunchao Zhang   /* create a vecscatter that just copies x to y */
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(x,isx,y,isy,&vscat));
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
11997929ea7SJunchao Zhang 
12097929ea7SJunchao Zhang   /* view y to check the result. y should be {0~63} */
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nBefore VecScatterRemap on SGToSG, SEQ vector y is:\n"));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
12397929ea7SJunchao Zhang 
12497929ea7SJunchao Zhang   /* now call the weird subroutine VecScatterRemap to slightly change the vecscatter.
12597929ea7SJunchao Zhang 
12697929ea7SJunchao Zhang      Create tomap as {32~63,0~31}. Originaly, we read from indices {0~64} of seq x to write to y. The remap
12797929ea7SJunchao Zhang      does indices[i] = tomap[indices[i]]. Therefore, after the remap, we read from indices{32~63,0~31} of seq x.
12897929ea7SJunchao Zhang    */
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&tomap));
13097929ea7SJunchao Zhang   for (i=0; i<n/2; i++) { tomap[i] = i+n/2; tomap[i+n/2] = i; };
1315f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterRemap(vscat,tomap,NULL));
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
13497929ea7SJunchao Zhang 
13597929ea7SJunchao Zhang   /* view y to check the result. y should be {32~63,0~31} */
1365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"After VecScatterRemap on SGToSG, SEQ vector y is:\n"));
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
13897929ea7SJunchao Zhang 
13997929ea7SJunchao Zhang   /* destroy everything before we recreate them in different types */
1405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(ix,iy));
1415f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1425f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isx));
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isy));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(tomap));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&vscat));
14797929ea7SJunchao Zhang 
14897929ea7SJunchao Zhang   /* ===================================================================================================
14997929ea7SJunchao Zhang      (3) test VecScatterRemap on a sequential general to sequential stride 1 (SGToSS_Stride1) vecscatter
15097929ea7SJunchao Zhang      ===================================================================================================
15197929ea7SJunchao Zhang    */
15297929ea7SJunchao Zhang   n = 64; /* long enough to trigger memcpy optimizations in local scatter */
15397929ea7SJunchao Zhang 
15497929ea7SJunchao Zhang   /* create two seq vectors x of length n, and y of length n/2 */
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&x));
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n/2,&y));
15797929ea7SJunchao Zhang 
15897929ea7SJunchao Zhang   /* Initialize x as {0~63} */
1595f80ce2aSJacob Faibussowitsch   for (i=0; i<n; i++) CHKERRQ(VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES));
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(x));
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(x));
16297929ea7SJunchao Zhang 
16397929ea7SJunchao Zhang   /* create a general index set isx = {0:63:2}, which actually is a stride IS with first=0, n=32, step=2,
16497929ea7SJunchao Zhang      but we use it as general and let PETSc detect the pattern and optimize it. */
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(n/2,&ix,n/2,&iy));
16697929ea7SJunchao Zhang   for (i=0; i<n/2; i++) ix[i] = i*2;
1675f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,n/2,ix,PETSC_COPY_VALUES,&isx));
16897929ea7SJunchao Zhang 
16997929ea7SJunchao Zhang   /* create a stride1 index set isy = {0~31}. We intentionally set the step to 1 to trigger optimizations */
1705f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateStride(PETSC_COMM_SELF,32,0,1,&isy));
17197929ea7SJunchao Zhang 
17297929ea7SJunchao Zhang   /* create a vecscatter that just copies even entries of x to y */
1735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(x,isx,y,isy,&vscat));
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
17697929ea7SJunchao Zhang 
17797929ea7SJunchao Zhang   /* view y to check the result. y should be {0:63:2} */
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nBefore VecScatterRemap on SGToSS_Stride1, SEQ vector y is:\n"));
1795f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
18097929ea7SJunchao Zhang 
18197929ea7SJunchao Zhang   /* now call the weird subroutine VecScatterRemap to slightly change the vecscatter.
18297929ea7SJunchao Zhang 
18397929ea7SJunchao Zhang      Create tomap as {32~63,0~31}. Originaly, we read from indices{0:63:2} of seq x to write to y. The remap
18497929ea7SJunchao Zhang      does indices[i] = tomap[indices[i]]. Therefore, after the remap, we read from indices{32:63:2,0:31:2} of seq x.
18597929ea7SJunchao Zhang    */
1865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&tomap));
18797929ea7SJunchao Zhang   for (i=0; i<n/2; i++) { tomap[i] = i+n/2; tomap[i+n/2] = i; };
1885f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterRemap(vscat,tomap,NULL));
1895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
1905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
19197929ea7SJunchao Zhang 
19297929ea7SJunchao Zhang   /* view y to check the result. y should be {32:63:2,0:31:2} */
1935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"After VecScatterRemap on SGToSS_Stride1, SEQ vector y is:\n"));
1945f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
19597929ea7SJunchao Zhang 
19697929ea7SJunchao Zhang   /* destroy everything before PetscFinalize */
1975f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(ix,iy));
1985f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
2005f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isx));
2015f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isy));
2025f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(tomap));
2035f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&vscat));
20497929ea7SJunchao Zhang 
205*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
206*b122ec5aSJacob Faibussowitsch   return 0;
20797929ea7SJunchao Zhang }
20897929ea7SJunchao Zhang 
20997929ea7SJunchao Zhang /*TEST
21097929ea7SJunchao Zhang 
21197929ea7SJunchao Zhang    test:
21297929ea7SJunchao Zhang       suffix: 1
21397929ea7SJunchao Zhang       nsize: 2
21497929ea7SJunchao Zhang       diff_args: -j
21597929ea7SJunchao Zhang       requires: double
21697929ea7SJunchao Zhang TEST*/
217