xref: /petsc/src/vec/is/sf/tests/ex15.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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   PetscErrorCode     ierr;
1497929ea7SJunchao Zhang   PetscInt           i,n,*ix,*iy,*tomap,start;
1597929ea7SJunchao Zhang   Vec                x,y;
1697929ea7SJunchao Zhang   PetscMPIInt        nproc,rank;
1797929ea7SJunchao Zhang   IS                 isx,isy;
1897929ea7SJunchao Zhang   const PetscInt     *ranges;
1997929ea7SJunchao Zhang   VecScatter         vscat;
2097929ea7SJunchao Zhang 
2197929ea7SJunchao Zhang   PetscFunctionBegin;
2297929ea7SJunchao Zhang   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
23ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&nproc);CHKERRMPI(ierr);
24ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
2597929ea7SJunchao Zhang 
26*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(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 */
3697929ea7SJunchao Zhang   ierr = VecCreateMPI(PETSC_COMM_WORLD,n,PETSC_DECIDE,&x);CHKERRQ(ierr);
3797929ea7SJunchao Zhang   ierr = VecDuplicate(x,&y);CHKERRQ(ierr);
3897929ea7SJunchao Zhang 
3997929ea7SJunchao Zhang   /* Initialize x as {0~127} */
4097929ea7SJunchao Zhang   ierr = VecGetOwnershipRanges(x,&ranges);CHKERRQ(ierr);
4197929ea7SJunchao Zhang   for (i=ranges[rank]; i<ranges[rank+1]; i++) { ierr = VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES);CHKERRQ(ierr); }
4297929ea7SJunchao Zhang   ierr = VecAssemblyBegin(x);CHKERRQ(ierr);
4397929ea7SJunchao Zhang   ierr = VecAssemblyEnd(x);CHKERRQ(ierr);
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    */
4997929ea7SJunchao Zhang   ierr = PetscMalloc2(n,&ix,n,&iy);CHKERRQ(ierr);
5097929ea7SJunchao Zhang   start = ranges[rank];
5197929ea7SJunchao Zhang   for (i=ranges[rank]; i<ranges[rank+1]; i++) ix[i-start] = i;
5297929ea7SJunchao Zhang   ierr = ISCreateGeneral(PETSC_COMM_WORLD,n,ix,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr);
5397929ea7SJunchao Zhang 
54dd400576SPatrick Sanan   if (rank == 0) { for (i=0; i<n; i++) iy[i] = i+32; }
5597929ea7SJunchao Zhang   else for (i=0; i<n/2; i++) { iy[i] = i+96; iy[i+n/2] = i; }
5697929ea7SJunchao Zhang 
5797929ea7SJunchao Zhang   ierr = ISCreateGeneral(PETSC_COMM_WORLD,n,iy,PETSC_COPY_VALUES,&isy);CHKERRQ(ierr);
5897929ea7SJunchao Zhang 
5997929ea7SJunchao Zhang   /* create a vecscatter that shifts x to the tail by quater periodically and puts the results in y */
6097929ea7SJunchao Zhang   ierr = VecScatterCreate(x,isx,y,isy,&vscat);CHKERRQ(ierr);
6197929ea7SJunchao Zhang   ierr = VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6297929ea7SJunchao Zhang   ierr = VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6397929ea7SJunchao Zhang 
6497929ea7SJunchao 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} */
6597929ea7SJunchao Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"Before VecScatterRemap on PtoP, MPI vector y is:\n");CHKERRQ(ierr);
6697929ea7SJunchao Zhang   ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
6797929ea7SJunchao Zhang 
6897929ea7SJunchao Zhang   /* now call the weird subroutine VecScatterRemap to slightly change the vecscatter. It changes where we read vector
6997929ea7SJunchao Zhang      x entries to send out, but does not change the communication pattern (i.e., send/recv pairs and msg lengths).
7097929ea7SJunchao Zhang 
7197929ea7SJunchao 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
7297929ea7SJunchao Zhang      does indices[i] = tomap[indices[i]]. Therefore, after the remap, we read from indices {32~63,0~31} of the local x.
7397929ea7SJunchao Zhang      isy is unchanged. So, we will shift x to {Q2,Q1,Q0,Q3}, that is {64~95,32~63,0~31,96~127}
7497929ea7SJunchao Zhang   */
7597929ea7SJunchao Zhang   ierr = PetscMalloc1(n,&tomap);CHKERRQ(ierr);
7697929ea7SJunchao Zhang   for (i=0; i<n/2; i++) { tomap[i] = i+n/2; tomap[i+n/2] = i; };
7797929ea7SJunchao Zhang   ierr = VecScatterRemap(vscat,tomap,NULL);CHKERRQ(ierr);
7897929ea7SJunchao Zhang   ierr = VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7997929ea7SJunchao Zhang   ierr = VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8097929ea7SJunchao Zhang 
8197929ea7SJunchao Zhang   /* view y to check the result. y should be {64~95,32~63,0~31,96~127} */
8297929ea7SJunchao Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"After VecScatterRemap on PtoP, MPI vector y is:\n");CHKERRQ(ierr);
8397929ea7SJunchao Zhang   ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
8497929ea7SJunchao Zhang 
8597929ea7SJunchao Zhang   /* destroy everything before we recreate them in different types */
8697929ea7SJunchao Zhang   ierr = PetscFree2(ix,iy);CHKERRQ(ierr);
8797929ea7SJunchao Zhang   ierr = VecDestroy(&x);CHKERRQ(ierr);
8897929ea7SJunchao Zhang   ierr = VecDestroy(&y);CHKERRQ(ierr);
8997929ea7SJunchao Zhang   ierr = ISDestroy(&isx);CHKERRQ(ierr);
9097929ea7SJunchao Zhang   ierr = ISDestroy(&isy);CHKERRQ(ierr);
9197929ea7SJunchao Zhang   ierr = PetscFree(tomap);CHKERRQ(ierr);
9297929ea7SJunchao Zhang   ierr = VecScatterDestroy(&vscat);CHKERRQ(ierr);
9397929ea7SJunchao Zhang 
9497929ea7SJunchao Zhang   /* ==========================================================================================
9597929ea7SJunchao Zhang      (2) test VecScatterRemap on a sequential general to sequential general (SGToSG) vecscatter
9697929ea7SJunchao Zhang      ==========================================================================================
9797929ea7SJunchao Zhang    */
9897929ea7SJunchao Zhang   n = 64; /* long enough to trigger memcpy optimizations in local scatter */
9997929ea7SJunchao Zhang 
10097929ea7SJunchao Zhang   /* create two seq vectors x, y of length n */
10197929ea7SJunchao Zhang   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&x);CHKERRQ(ierr);
10297929ea7SJunchao Zhang   ierr = VecDuplicate(x,&y);CHKERRQ(ierr);
10397929ea7SJunchao Zhang 
10497929ea7SJunchao Zhang   /* Initialize x as {0~63} */
10597929ea7SJunchao Zhang   for (i=0; i<n; i++) { ierr = VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES);CHKERRQ(ierr); }
10697929ea7SJunchao Zhang   ierr = VecAssemblyBegin(x);CHKERRQ(ierr);
10797929ea7SJunchao Zhang   ierr = VecAssemblyEnd(x);CHKERRQ(ierr);
10897929ea7SJunchao Zhang 
10997929ea7SJunchao Zhang   /* create two general index sets isx = isy = {0~63}, which are sequential, but we use them as
11097929ea7SJunchao Zhang      general and let PETSc detect the pattern and optimize it */
11197929ea7SJunchao Zhang   ierr = PetscMalloc2(n,&ix,n,&iy);CHKERRQ(ierr);
11297929ea7SJunchao Zhang   for (i=0; i<n; i++) ix[i] = i;
11397929ea7SJunchao Zhang   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ix,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr);
11497929ea7SJunchao Zhang   ierr = ISDuplicate(isx,&isy);CHKERRQ(ierr);
11597929ea7SJunchao Zhang 
11697929ea7SJunchao Zhang   /* create a vecscatter that just copies x to y */
11797929ea7SJunchao Zhang   ierr = VecScatterCreate(x,isx,y,isy,&vscat);CHKERRQ(ierr);
11897929ea7SJunchao Zhang   ierr = VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11997929ea7SJunchao Zhang   ierr = VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12097929ea7SJunchao Zhang 
12197929ea7SJunchao Zhang   /* view y to check the result. y should be {0~63} */
12297929ea7SJunchao Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"\nBefore VecScatterRemap on SGToSG, SEQ vector y is:\n");CHKERRQ(ierr);
12397929ea7SJunchao Zhang   ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
12497929ea7SJunchao Zhang 
12597929ea7SJunchao Zhang   /* now call the weird subroutine VecScatterRemap to slightly change the vecscatter.
12697929ea7SJunchao Zhang 
12797929ea7SJunchao Zhang      Create tomap as {32~63,0~31}. Originaly, we read from indices {0~64} of seq x to write to y. The remap
12897929ea7SJunchao Zhang      does indices[i] = tomap[indices[i]]. Therefore, after the remap, we read from indices{32~63,0~31} of seq x.
12997929ea7SJunchao Zhang    */
13097929ea7SJunchao Zhang   ierr = PetscMalloc1(n,&tomap);CHKERRQ(ierr);
13197929ea7SJunchao Zhang   for (i=0; i<n/2; i++) { tomap[i] = i+n/2; tomap[i+n/2] = i; };
13297929ea7SJunchao Zhang   ierr = VecScatterRemap(vscat,tomap,NULL);CHKERRQ(ierr);
13397929ea7SJunchao Zhang   ierr = VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13497929ea7SJunchao Zhang   ierr = VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13597929ea7SJunchao Zhang 
13697929ea7SJunchao Zhang   /* view y to check the result. y should be {32~63,0~31} */
13797929ea7SJunchao Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"After VecScatterRemap on SGToSG, SEQ vector y is:\n");CHKERRQ(ierr);
13897929ea7SJunchao Zhang   ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
13997929ea7SJunchao Zhang 
14097929ea7SJunchao Zhang   /* destroy everything before we recreate them in different types */
14197929ea7SJunchao Zhang   ierr = PetscFree2(ix,iy);CHKERRQ(ierr);
14297929ea7SJunchao Zhang   ierr = VecDestroy(&x);CHKERRQ(ierr);
14397929ea7SJunchao Zhang   ierr = VecDestroy(&y);CHKERRQ(ierr);
14497929ea7SJunchao Zhang   ierr = ISDestroy(&isx);CHKERRQ(ierr);
14597929ea7SJunchao Zhang   ierr = ISDestroy(&isy);CHKERRQ(ierr);
14697929ea7SJunchao Zhang   ierr = PetscFree(tomap);CHKERRQ(ierr);
14797929ea7SJunchao Zhang   ierr = VecScatterDestroy(&vscat);CHKERRQ(ierr);
14897929ea7SJunchao Zhang 
14997929ea7SJunchao Zhang   /* ===================================================================================================
15097929ea7SJunchao Zhang      (3) test VecScatterRemap on a sequential general to sequential stride 1 (SGToSS_Stride1) vecscatter
15197929ea7SJunchao Zhang      ===================================================================================================
15297929ea7SJunchao Zhang    */
15397929ea7SJunchao Zhang   n = 64; /* long enough to trigger memcpy optimizations in local scatter */
15497929ea7SJunchao Zhang 
15597929ea7SJunchao Zhang   /* create two seq vectors x of length n, and y of length n/2 */
15697929ea7SJunchao Zhang   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&x);CHKERRQ(ierr);
15797929ea7SJunchao Zhang   ierr = VecCreateSeq(PETSC_COMM_SELF,n/2,&y);CHKERRQ(ierr);
15897929ea7SJunchao Zhang 
15997929ea7SJunchao Zhang   /* Initialize x as {0~63} */
16097929ea7SJunchao Zhang   for (i=0; i<n; i++) { ierr = VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES);CHKERRQ(ierr); }
16197929ea7SJunchao Zhang   ierr = VecAssemblyBegin(x);CHKERRQ(ierr);
16297929ea7SJunchao Zhang   ierr = VecAssemblyEnd(x);CHKERRQ(ierr);
16397929ea7SJunchao Zhang 
16497929ea7SJunchao Zhang   /* create a general index set isx = {0:63:2}, which actually is a stride IS with first=0, n=32, step=2,
16597929ea7SJunchao Zhang      but we use it as general and let PETSc detect the pattern and optimize it. */
16697929ea7SJunchao Zhang   ierr = PetscMalloc2(n/2,&ix,n/2,&iy);CHKERRQ(ierr);
16797929ea7SJunchao Zhang   for (i=0; i<n/2; i++) ix[i] = i*2;
16897929ea7SJunchao Zhang   ierr = ISCreateGeneral(PETSC_COMM_SELF,n/2,ix,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr);
16997929ea7SJunchao Zhang 
17097929ea7SJunchao Zhang   /* create a stride1 index set isy = {0~31}. We intentionally set the step to 1 to trigger optimizations */
17197929ea7SJunchao Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,32,0,1,&isy);CHKERRQ(ierr);
17297929ea7SJunchao Zhang 
17397929ea7SJunchao Zhang   /* create a vecscatter that just copies even entries of x to y */
17497929ea7SJunchao Zhang   ierr = VecScatterCreate(x,isx,y,isy,&vscat);CHKERRQ(ierr);
17597929ea7SJunchao Zhang   ierr = VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
17697929ea7SJunchao Zhang   ierr = VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
17797929ea7SJunchao Zhang 
17897929ea7SJunchao Zhang   /* view y to check the result. y should be {0:63:2} */
17997929ea7SJunchao Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"\nBefore VecScatterRemap on SGToSS_Stride1, SEQ vector y is:\n");CHKERRQ(ierr);
18097929ea7SJunchao Zhang   ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
18197929ea7SJunchao Zhang 
18297929ea7SJunchao Zhang   /* now call the weird subroutine VecScatterRemap to slightly change the vecscatter.
18397929ea7SJunchao Zhang 
18497929ea7SJunchao 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
18597929ea7SJunchao Zhang      does indices[i] = tomap[indices[i]]. Therefore, after the remap, we read from indices{32:63:2,0:31:2} of seq x.
18697929ea7SJunchao Zhang    */
18797929ea7SJunchao Zhang   ierr = PetscMalloc1(n,&tomap);CHKERRQ(ierr);
18897929ea7SJunchao Zhang   for (i=0; i<n/2; i++) { tomap[i] = i+n/2; tomap[i+n/2] = i; };
18997929ea7SJunchao Zhang   ierr = VecScatterRemap(vscat,tomap,NULL);CHKERRQ(ierr);
19097929ea7SJunchao Zhang   ierr = VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
19197929ea7SJunchao Zhang   ierr = VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
19297929ea7SJunchao Zhang 
19397929ea7SJunchao Zhang   /* view y to check the result. y should be {32:63:2,0:31:2} */
19497929ea7SJunchao Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"After VecScatterRemap on SGToSS_Stride1, SEQ vector y is:\n");CHKERRQ(ierr);
19597929ea7SJunchao Zhang   ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
19697929ea7SJunchao Zhang 
19797929ea7SJunchao Zhang   /* destroy everything before PetscFinalize */
19897929ea7SJunchao Zhang   ierr = PetscFree2(ix,iy);CHKERRQ(ierr);
19997929ea7SJunchao Zhang   ierr = VecDestroy(&x);CHKERRQ(ierr);
20097929ea7SJunchao Zhang   ierr = VecDestroy(&y);CHKERRQ(ierr);
20197929ea7SJunchao Zhang   ierr = ISDestroy(&isx);CHKERRQ(ierr);
20297929ea7SJunchao Zhang   ierr = ISDestroy(&isy);CHKERRQ(ierr);
20397929ea7SJunchao Zhang   ierr = PetscFree(tomap);CHKERRQ(ierr);
20497929ea7SJunchao Zhang   ierr = VecScatterDestroy(&vscat);CHKERRQ(ierr);
20597929ea7SJunchao Zhang 
20697929ea7SJunchao Zhang   ierr = PetscFinalize();
20797929ea7SJunchao Zhang   return ierr;
20897929ea7SJunchao Zhang }
20997929ea7SJunchao Zhang 
21097929ea7SJunchao Zhang /*TEST
21197929ea7SJunchao Zhang 
21297929ea7SJunchao Zhang    test:
21397929ea7SJunchao Zhang       suffix: 1
21497929ea7SJunchao Zhang       nsize: 2
21597929ea7SJunchao Zhang       diff_args: -j
21697929ea7SJunchao Zhang       requires: double
21797929ea7SJunchao Zhang TEST*/
21897929ea7SJunchao Zhang 
219