xref: /petsc/src/vec/is/sf/tests/ex9.c (revision ffc4695bcb29f4b022f59a5fd6bc99fc280ff6d8)
197929ea7SJunchao Zhang static char help[]= "This example shows 1) how to transfer vectors from a parent communicator to vectors on a child communicator and vice versa;\n\
297929ea7SJunchao Zhang   2) how to transfer vectors from a subcommunicator to vectors on another subcommunicator. The two subcommunicators are not\n\
397929ea7SJunchao Zhang   required to cover all processes in PETSC_COMM_WORLD; 3) how to copy a vector from a parent communicator to vectors on its child communicators.\n\n";
497929ea7SJunchao Zhang 
597929ea7SJunchao Zhang #include <petscvec.h>
697929ea7SJunchao Zhang int main(int argc,char **argv)
797929ea7SJunchao Zhang {
897929ea7SJunchao Zhang   PetscErrorCode ierr;
997929ea7SJunchao Zhang   PetscMPIInt    nproc,grank,mycolor;
1097929ea7SJunchao Zhang   PetscInt       i,n,N = 20,low,high;
1197929ea7SJunchao Zhang   MPI_Comm       subcomm;
1297929ea7SJunchao Zhang   Vec            x,yg; /* global vectors on PETSC_COMM_WORLD */
1397929ea7SJunchao Zhang   VecScatter     vscat;
1497929ea7SJunchao Zhang   IS             ix,iy;
1597929ea7SJunchao Zhang   PetscBool      world2sub  = PETSC_FALSE;  /* Copy a vector from WORLD to a subcomm? */
1697929ea7SJunchao Zhang   PetscBool      sub2sub    = PETSC_FALSE;  /* Copy a vector from a subcomm to another subcomm? */
1797929ea7SJunchao Zhang   PetscBool      world2subs = PETSC_FALSE;  /* Copy a vector from WORLD to multiple subcomms? */
1897929ea7SJunchao Zhang 
1997929ea7SJunchao Zhang   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
20*ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&nproc);CHKERRMPI(ierr);
21*ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&grank);CHKERRMPI(ierr);
2297929ea7SJunchao Zhang 
2397929ea7SJunchao Zhang   if (nproc < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_SIZ,"This test must have at least two processes to run");
2497929ea7SJunchao Zhang 
2597929ea7SJunchao Zhang   ierr = PetscOptionsGetBool(NULL,0,"-world2sub",&world2sub,NULL);CHKERRQ(ierr);
2697929ea7SJunchao Zhang   ierr = PetscOptionsGetBool(NULL,0,"-sub2sub",&sub2sub,NULL);CHKERRQ(ierr);
2797929ea7SJunchao Zhang   ierr = PetscOptionsGetBool(NULL,0,"-world2subs",&world2subs,NULL);CHKERRQ(ierr);
2897929ea7SJunchao Zhang 
2997929ea7SJunchao Zhang   /* Split PETSC_COMM_WORLD into three subcomms. Each process can only see the subcomm it belongs to */
3097929ea7SJunchao Zhang   mycolor = grank % 3;
3197929ea7SJunchao Zhang   ierr    = MPI_Comm_split(PETSC_COMM_WORLD,mycolor,grank,&subcomm);CHKERRQ(ierr);
3297929ea7SJunchao Zhang 
3397929ea7SJunchao Zhang   /*===========================================================================
3497929ea7SJunchao Zhang    *  Transfer a vector x defined on PETSC_COMM_WORLD to a vector y defined on
3597929ea7SJunchao Zhang    *  a subcommunicator of PETSC_COMM_WORLD and vice versa.
3697929ea7SJunchao Zhang    *===========================================================================*/
3797929ea7SJunchao Zhang   if (world2sub) {
3897929ea7SJunchao Zhang     ierr = VecCreateMPI(PETSC_COMM_WORLD,PETSC_DECIDE,N,&x);CHKERRQ(ierr);
3997929ea7SJunchao Zhang     ierr = PetscObjectSetName((PetscObject)x,"x_commworld");CHKERRQ(ierr); /* Give a name to view x clearly */
4097929ea7SJunchao Zhang 
4197929ea7SJunchao Zhang     /* Initialize x to [-0.0, -1.0, -2.0, ..., -19.0] */
4297929ea7SJunchao Zhang     ierr = VecGetOwnershipRange(x,&low,&high);CHKERRQ(ierr);
4397929ea7SJunchao Zhang     for (i=low; i<high; i++) {
4497929ea7SJunchao Zhang       PetscScalar val = -i;
4597929ea7SJunchao Zhang       ierr = VecSetValue(x,i,val,INSERT_VALUES);CHKERRQ(ierr);
4697929ea7SJunchao Zhang     }
4797929ea7SJunchao Zhang     ierr = VecAssemblyBegin(x);CHKERRQ(ierr);
4897929ea7SJunchao Zhang     ierr = VecAssemblyEnd(x);CHKERRQ(ierr);
4997929ea7SJunchao Zhang 
5097929ea7SJunchao Zhang     /* Transfer x to a vector y only defined on subcomm0 and vice versa */
5197929ea7SJunchao Zhang     if (mycolor == 0) { /* subcomm0 contains ranks 0, 3, 6, ... in PETSC_COMM_WORLD */
5297929ea7SJunchao Zhang       Vec         y;
5397929ea7SJunchao Zhang       PetscScalar *yvalue;
5497929ea7SJunchao Zhang 
5597929ea7SJunchao Zhang       ierr = VecCreateMPI(subcomm,PETSC_DECIDE,N,&y);CHKERRQ(ierr);
5697929ea7SJunchao Zhang       ierr = PetscObjectSetName((PetscObject)y,"y_subcomm_0");CHKERRQ(ierr); /* Give a name to view y clearly */
5797929ea7SJunchao Zhang       ierr = VecGetLocalSize(y,&n);CHKERRQ(ierr);
5897929ea7SJunchao Zhang       ierr = VecGetArray(y,&yvalue);CHKERRQ(ierr);
5997929ea7SJunchao Zhang       /* Create yg on PETSC_COMM_WORLD and alias yg with y. They share the memory pointed by yvalue.
6097929ea7SJunchao Zhang         Note this is a collective call. All processes have to call it and supply consistent N.
6197929ea7SJunchao Zhang       */
6297929ea7SJunchao Zhang       ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,N,yvalue,&yg);CHKERRQ(ierr);
6397929ea7SJunchao Zhang 
6497929ea7SJunchao Zhang       /* Create an identity map that makes yg[i] = x[i], i=0..N-1 */
6597929ea7SJunchao Zhang       ierr = VecGetOwnershipRange(yg,&low,&high);CHKERRQ(ierr); /* low, high are global indices */
6697929ea7SJunchao Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,high-low,low,1,&ix);CHKERRQ(ierr);
6797929ea7SJunchao Zhang       ierr = ISDuplicate(ix,&iy);CHKERRQ(ierr);
6897929ea7SJunchao Zhang 
6997929ea7SJunchao Zhang       /* Union of ix's on subcomm0 covers the full range of [0,N) */
7097929ea7SJunchao Zhang       ierr = VecScatterCreate(x,ix,yg,iy,&vscat);CHKERRQ(ierr);
7197929ea7SJunchao Zhang       ierr = VecScatterBegin(vscat,x,yg,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7297929ea7SJunchao Zhang       ierr = VecScatterEnd(vscat,x,yg,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7397929ea7SJunchao Zhang 
7497929ea7SJunchao Zhang       /* Once yg got the data from x, we return yvalue to y so that we can use y in other operations.
7597929ea7SJunchao Zhang         VecGetArray must be paired with VecRestoreArray.
7697929ea7SJunchao Zhang       */
7797929ea7SJunchao Zhang       ierr = VecRestoreArray(y,&yvalue);CHKERRQ(ierr);
7897929ea7SJunchao Zhang 
7997929ea7SJunchao Zhang       /* Libraries on subcomm0 can safely use y now, for example, view and scale it */
8097929ea7SJunchao Zhang       ierr = VecView(y,PETSC_VIEWER_STDOUT_(subcomm));CHKERRQ(ierr);
8197929ea7SJunchao Zhang       ierr = VecScale(y,2.0);CHKERRQ(ierr);
8297929ea7SJunchao Zhang 
8397929ea7SJunchao Zhang       /* Send the new y back to x */
8497929ea7SJunchao Zhang       ierr = VecGetArray(y,&yvalue);CHKERRQ(ierr); /* If VecScale is done on GPU, Petsc will prepare a valid yvalue for access */
8597929ea7SJunchao Zhang       /* Supply new yvalue to yg without memory copying */
8697929ea7SJunchao Zhang       ierr = VecPlaceArray(yg,yvalue);CHKERRQ(ierr);
8797929ea7SJunchao Zhang       ierr = VecScatterBegin(vscat,yg,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8897929ea7SJunchao Zhang       ierr = VecScatterEnd(vscat,yg,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8997929ea7SJunchao Zhang       ierr = VecResetArray(yg);CHKERRQ(ierr);
9097929ea7SJunchao Zhang       ierr = VecRestoreArray(y,&yvalue);CHKERRQ(ierr);
9197929ea7SJunchao Zhang 
9297929ea7SJunchao Zhang       ierr = VecDestroy(&y);CHKERRQ(ierr);
9397929ea7SJunchao Zhang     } else {
9497929ea7SJunchao Zhang       /* Ranks outside of subcomm0 do not supply values to yg */
9597929ea7SJunchao Zhang       ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,0/*n*/,N,NULL,&yg);CHKERRQ(ierr);
9697929ea7SJunchao Zhang 
9797929ea7SJunchao Zhang       /* Ranks in subcomm0 already specified the full range of the identity map. The remaining
9897929ea7SJunchao Zhang         ranks just need to create empty ISes to cheat VecScatterCreate.
9997929ea7SJunchao Zhang       */
10097929ea7SJunchao Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&ix);CHKERRQ(ierr);
10197929ea7SJunchao Zhang       ierr = ISDuplicate(ix,&iy);CHKERRQ(ierr);
10297929ea7SJunchao Zhang 
10397929ea7SJunchao Zhang       ierr = VecScatterCreate(x,ix,yg,iy,&vscat);CHKERRQ(ierr);
10497929ea7SJunchao Zhang       ierr = VecScatterBegin(vscat,x,yg,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10597929ea7SJunchao Zhang       ierr = VecScatterEnd(vscat,x,yg,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10697929ea7SJunchao Zhang 
10797929ea7SJunchao Zhang       /* Send the new y back to x. Ranks outside of subcomm0 actually have nothing to send.
10897929ea7SJunchao Zhang         But they have to call VecScatterBegin/End since these routines are collective.
10997929ea7SJunchao Zhang       */
11097929ea7SJunchao Zhang       ierr = VecScatterBegin(vscat,yg,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11197929ea7SJunchao Zhang       ierr = VecScatterEnd(vscat,yg,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11297929ea7SJunchao Zhang     }
11397929ea7SJunchao Zhang 
11497929ea7SJunchao Zhang     ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
11597929ea7SJunchao Zhang     ierr = ISDestroy(&ix);CHKERRQ(ierr);
11697929ea7SJunchao Zhang     ierr = ISDestroy(&iy);CHKERRQ(ierr);
11797929ea7SJunchao Zhang     ierr = VecDestroy(&x);CHKERRQ(ierr);
11897929ea7SJunchao Zhang     ierr = VecDestroy(&yg);CHKERRQ(ierr);
11997929ea7SJunchao Zhang     ierr = VecScatterDestroy(&vscat);CHKERRQ(ierr);
12097929ea7SJunchao Zhang   } /* world2sub */
12197929ea7SJunchao Zhang 
12297929ea7SJunchao Zhang   /*===========================================================================
12397929ea7SJunchao Zhang    *  Transfer a vector x defined on subcomm0 to a vector y defined on
12497929ea7SJunchao Zhang    *  subcomm1. The two subcomms are not overlapping and their union is
12597929ea7SJunchao Zhang    *  not necessarily equal to PETSC_COMM_WORLD.
12697929ea7SJunchao Zhang    *===========================================================================*/
12797929ea7SJunchao Zhang   if (sub2sub) {
12897929ea7SJunchao Zhang     if (mycolor == 0) {
12997929ea7SJunchao Zhang       /* Intentionally declare N as a local variable so that processes in subcomm1 do not know its value */
13097929ea7SJunchao Zhang       PetscInt          n,N = 22;
13197929ea7SJunchao Zhang       Vec               x,xg,yg;
13297929ea7SJunchao Zhang       IS                ix,iy;
13397929ea7SJunchao Zhang       VecScatter        vscat;
13497929ea7SJunchao Zhang       const PetscScalar *xvalue;
13597929ea7SJunchao Zhang       MPI_Comm          intercomm,parentcomm;
13697929ea7SJunchao Zhang       PetscMPIInt       lrank;
13797929ea7SJunchao Zhang 
138*ffc4695bSBarry Smith       ierr = MPI_Comm_rank(subcomm,&lrank);CHKERRMPI(ierr);
13997929ea7SJunchao Zhang       ierr = VecCreateMPI(subcomm,PETSC_DECIDE,N,&x);CHKERRQ(ierr); /* x is on subcomm */
14097929ea7SJunchao Zhang       ierr = VecGetOwnershipRange(x,&low,&high);CHKERRQ(ierr);
14197929ea7SJunchao Zhang 
14297929ea7SJunchao Zhang       /* initialize x = [0.0, 1.0, 2.0, ..., 21.0] */
14397929ea7SJunchao Zhang       for (i=low; i<high; i++) {
14497929ea7SJunchao Zhang         PetscScalar val = i;
14597929ea7SJunchao Zhang         ierr = VecSetValue(x,i,val,INSERT_VALUES);CHKERRQ(ierr);
14697929ea7SJunchao Zhang       }
14797929ea7SJunchao Zhang       ierr = VecAssemblyBegin(x);CHKERRQ(ierr);
14897929ea7SJunchao Zhang       ierr = VecAssemblyEnd(x);CHKERRQ(ierr);
14997929ea7SJunchao Zhang 
150*ffc4695bSBarry Smith       ierr = MPI_Intercomm_create(subcomm,0,PETSC_COMM_WORLD/*peer_comm*/,1,100/*tag*/,&intercomm);CHKERRMPI(ierr);
15197929ea7SJunchao Zhang 
15297929ea7SJunchao Zhang       /* Tell rank 0 of subcomm1 the global size of x */
153*ffc4695bSBarry Smith       if (!lrank) {ierr = MPI_Send(&N,1,MPIU_INT,0/*receiver's rank in remote comm, i.e., subcomm1*/,200/*tag*/,intercomm);CHKERRMPI(ierr);}
15497929ea7SJunchao Zhang 
15597929ea7SJunchao Zhang       /* Create an intracomm Petsc can work on. Ranks in subcomm0 are ordered before ranks in subcomm1 in parentcomm.
15697929ea7SJunchao Zhang         But this order actually does not matter, since what we care is vector y, which is defined on subcomm1.
15797929ea7SJunchao Zhang       */
158*ffc4695bSBarry Smith       ierr = MPI_Intercomm_merge(intercomm,0/*low*/,&parentcomm);CHKERRMPI(ierr);
15997929ea7SJunchao Zhang 
16097929ea7SJunchao Zhang       /* Create a vector xg on parentcomm, which shares memory with x */
16197929ea7SJunchao Zhang       ierr = VecGetArrayRead(x,&xvalue);CHKERRQ(ierr);
16297929ea7SJunchao Zhang       ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr);
16397929ea7SJunchao Zhang       ierr = VecCreateMPIWithArray(parentcomm,1,n,N,xvalue,&xg);CHKERRQ(ierr);
16497929ea7SJunchao Zhang 
16597929ea7SJunchao Zhang       /* Ranks in subcomm 0 have nothing on yg, so they simply have n=0, array=NULL */
16697929ea7SJunchao Zhang       ierr = VecCreateMPIWithArray(parentcomm,1,0/*n*/,N,NULL/*array*/,&yg);CHKERRQ(ierr);
16797929ea7SJunchao Zhang 
16897929ea7SJunchao Zhang       /* Create the vecscatter, which does identity map by setting yg[i] = xg[i], i=0..N-1. */
16997929ea7SJunchao Zhang       ierr = VecGetOwnershipRange(xg,&low,&high);CHKERRQ(ierr); /* low, high are global indices of xg */
17097929ea7SJunchao Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,high-low,low,1,&ix);CHKERRQ(ierr);
17197929ea7SJunchao Zhang       ierr = ISDuplicate(ix,&iy);CHKERRQ(ierr);
17297929ea7SJunchao Zhang       ierr = VecScatterCreate(xg,ix,yg,iy,&vscat);CHKERRQ(ierr);
17397929ea7SJunchao Zhang 
17497929ea7SJunchao Zhang       /* Scatter values from xg to yg */
17597929ea7SJunchao Zhang       ierr = VecScatterBegin(vscat,xg,yg,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
17697929ea7SJunchao Zhang       ierr = VecScatterEnd(vscat,xg,yg,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
17797929ea7SJunchao Zhang 
17897929ea7SJunchao Zhang       /* After the VecScatter is done, xg is idle so we can safely return xvalue to x */
17997929ea7SJunchao Zhang       ierr = VecRestoreArrayRead(x,&xvalue);CHKERRQ(ierr);
18097929ea7SJunchao Zhang       ierr = VecDestroy(&x);CHKERRQ(ierr);
18197929ea7SJunchao Zhang       ierr = ISDestroy(&ix);CHKERRQ(ierr);
18297929ea7SJunchao Zhang       ierr = ISDestroy(&iy);CHKERRQ(ierr);
18397929ea7SJunchao Zhang       ierr = VecDestroy(&xg);CHKERRQ(ierr);
18497929ea7SJunchao Zhang       ierr = VecDestroy(&yg);CHKERRQ(ierr);
18597929ea7SJunchao Zhang       ierr = VecScatterDestroy(&vscat);CHKERRQ(ierr);
186*ffc4695bSBarry Smith       ierr = MPI_Comm_free(&intercomm);CHKERRMPI(ierr);
187*ffc4695bSBarry Smith       ierr = MPI_Comm_free(&parentcomm);CHKERRMPI(ierr);
18897929ea7SJunchao Zhang     } else if (mycolor == 1) { /* subcomm 1, containing ranks 1, 4, 7, ... in PETSC_COMM_WORLD */
18997929ea7SJunchao Zhang       PetscInt    n,N;
19097929ea7SJunchao Zhang       Vec         y,xg,yg;
19197929ea7SJunchao Zhang       IS          ix,iy;
19297929ea7SJunchao Zhang       VecScatter  vscat;
19397929ea7SJunchao Zhang       PetscScalar *yvalue;
19497929ea7SJunchao Zhang       MPI_Comm    intercomm,parentcomm;
19597929ea7SJunchao Zhang       PetscMPIInt lrank;
19697929ea7SJunchao Zhang 
197*ffc4695bSBarry Smith       ierr = MPI_Comm_rank(subcomm,&lrank);CHKERRMPI(ierr);
198*ffc4695bSBarry Smith       ierr = MPI_Intercomm_create(subcomm,0,PETSC_COMM_WORLD/*peer_comm*/,0/*remote_leader*/,100/*tag*/,&intercomm);CHKERRMPI(ierr);
19997929ea7SJunchao Zhang 
20097929ea7SJunchao Zhang       /* Two rank-0 are talking */
201*ffc4695bSBarry Smith       if (!lrank) {ierr = MPI_Recv(&N,1,MPIU_INT,0/*sender's rank in remote comm, i.e. subcomm0*/,200/*tag*/,intercomm,MPI_STATUS_IGNORE);CHKERRMPI(ierr);}
20297929ea7SJunchao Zhang       /* Rank 0 of subcomm1 bcasts N to its members */
203*ffc4695bSBarry Smith       ierr = MPI_Bcast(&N,1,MPIU_INT,0/*local root*/,subcomm);CHKERRMPI(ierr);
20497929ea7SJunchao Zhang 
20597929ea7SJunchao Zhang       /* Create a intracomm Petsc can work on */
206*ffc4695bSBarry Smith       ierr = MPI_Intercomm_merge(intercomm,1/*high*/,&parentcomm);CHKERRMPI(ierr);
20797929ea7SJunchao Zhang 
20897929ea7SJunchao Zhang       /* Ranks in subcomm1 have nothing on xg, so they simply have n=0, array=NULL.*/
20997929ea7SJunchao Zhang       ierr = VecCreateMPIWithArray(parentcomm,1/*bs*/,0/*n*/,N,NULL/*array*/,&xg);CHKERRQ(ierr);
21097929ea7SJunchao Zhang 
21197929ea7SJunchao Zhang       ierr = VecCreateMPI(subcomm,PETSC_DECIDE,N,&y);CHKERRQ(ierr);
21297929ea7SJunchao Zhang       ierr = PetscObjectSetName((PetscObject)y,"y_subcomm_1");CHKERRQ(ierr); /* Give a name to view y clearly */
21397929ea7SJunchao Zhang       ierr = VecGetLocalSize(y,&n);CHKERRQ(ierr);
21497929ea7SJunchao Zhang       ierr = VecGetArray(y,&yvalue);CHKERRQ(ierr);
21597929ea7SJunchao Zhang       /* Create a vector yg on parentcomm, which shares memory with y. xg and yg must be
21697929ea7SJunchao Zhang         created in the same order in subcomm0/1. For example, we can not reverse the order of
21797929ea7SJunchao Zhang         creating xg and yg in subcomm1.
21897929ea7SJunchao Zhang       */
21997929ea7SJunchao Zhang       ierr = VecCreateMPIWithArray(parentcomm,1/*bs*/,n,N,yvalue,&yg);CHKERRQ(ierr);
22097929ea7SJunchao Zhang 
22197929ea7SJunchao Zhang       /* Ranks in subcomm0 already specified the full range of the identity map.
22297929ea7SJunchao Zhang         ranks in subcomm1 just need to create empty ISes to cheat VecScatterCreate.
22397929ea7SJunchao Zhang       */
22497929ea7SJunchao Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&ix);CHKERRQ(ierr);
22597929ea7SJunchao Zhang       ierr = ISDuplicate(ix,&iy);CHKERRQ(ierr);
22697929ea7SJunchao Zhang       ierr = VecScatterCreate(xg,ix,yg,iy,&vscat);CHKERRQ(ierr);
22797929ea7SJunchao Zhang 
22897929ea7SJunchao Zhang       /* Scatter values from xg to yg */
22997929ea7SJunchao Zhang       ierr = VecScatterBegin(vscat,xg,yg,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
23097929ea7SJunchao Zhang       ierr = VecScatterEnd(vscat,xg,yg,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
23197929ea7SJunchao Zhang 
23297929ea7SJunchao Zhang       /* After the VecScatter is done, values in yg are available. y is our interest, so we return yvalue to y */
23397929ea7SJunchao Zhang       ierr = VecRestoreArray(y,&yvalue);CHKERRQ(ierr);
23497929ea7SJunchao Zhang 
23597929ea7SJunchao Zhang       /* Libraries on subcomm1 can safely use y now, for example, view it */
23697929ea7SJunchao Zhang       ierr = VecView(y,PETSC_VIEWER_STDOUT_(subcomm));CHKERRQ(ierr);
23797929ea7SJunchao Zhang 
23897929ea7SJunchao Zhang       ierr = VecDestroy(&y);CHKERRQ(ierr);
23997929ea7SJunchao Zhang       ierr = ISDestroy(&ix);CHKERRQ(ierr);
24097929ea7SJunchao Zhang       ierr = ISDestroy(&iy);CHKERRQ(ierr);
24197929ea7SJunchao Zhang       ierr = VecDestroy(&xg);CHKERRQ(ierr);
24297929ea7SJunchao Zhang       ierr = VecDestroy(&yg);CHKERRQ(ierr);
24397929ea7SJunchao Zhang       ierr = VecScatterDestroy(&vscat);CHKERRQ(ierr);
244*ffc4695bSBarry Smith       ierr = MPI_Comm_free(&intercomm);CHKERRMPI(ierr);
245*ffc4695bSBarry Smith       ierr = MPI_Comm_free(&parentcomm);CHKERRMPI(ierr);
24697929ea7SJunchao Zhang     } else if (mycolor == 2) { /* subcomm2 */
24797929ea7SJunchao Zhang       /* Processes in subcomm2 do not participate in the VecScatter. They can freely do unrelated things on subcomm2 */
24897929ea7SJunchao Zhang     }
24997929ea7SJunchao Zhang   } /* sub2sub */
25097929ea7SJunchao Zhang 
25197929ea7SJunchao Zhang   /*===========================================================================
25297929ea7SJunchao Zhang    *  Copy a vector x defined on PETSC_COMM_WORLD to vectors y defined on
25397929ea7SJunchao Zhang    *  every subcommunicator of PETSC_COMM_WORLD. We could use multiple transfers
25497929ea7SJunchao Zhang    *  as we did in case 1, but that is not efficient. Instead, we use one vecscatter
25597929ea7SJunchao Zhang    *  to achieve that.
25697929ea7SJunchao Zhang    *===========================================================================*/
25797929ea7SJunchao Zhang   if (world2subs) {
25897929ea7SJunchao Zhang     Vec         y,yg;
25997929ea7SJunchao Zhang     PetscInt    n,N=15,xstart,ystart,low,high;
26097929ea7SJunchao Zhang     PetscScalar *yvalue;
26197929ea7SJunchao Zhang 
26297929ea7SJunchao Zhang     /* Initialize x to [0, 1, 2, 3, ..., N-1] */
26397929ea7SJunchao Zhang     ierr = VecCreateMPI(PETSC_COMM_WORLD,PETSC_DECIDE,N,&x);CHKERRQ(ierr);
26497929ea7SJunchao Zhang     ierr = VecGetOwnershipRange(x,&low,&high);CHKERRQ(ierr);
26597929ea7SJunchao Zhang     for (i=low; i<high; i++) {ierr = VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES);CHKERRQ(ierr);}
26697929ea7SJunchao Zhang     ierr = VecAssemblyBegin(x);CHKERRQ(ierr);
26797929ea7SJunchao Zhang     ierr = VecAssemblyEnd(x);CHKERRQ(ierr);
26897929ea7SJunchao Zhang 
26997929ea7SJunchao Zhang     /* Every subcomm has a y as long as x */
27097929ea7SJunchao Zhang     ierr = VecCreateMPI(subcomm,PETSC_DECIDE,N,&y);CHKERRQ(ierr);
27197929ea7SJunchao Zhang     ierr = VecGetLocalSize(y,&n);CHKERRQ(ierr);
27297929ea7SJunchao Zhang 
27397929ea7SJunchao Zhang     /* Create a global vector yg on PETSC_COMM_WORLD using y's memory. yg's global size = N*(number of subcommunicators).
27497929ea7SJunchao Zhang        Eeach rank in subcomms contributes a piece to construct the global yg. Keep in mind that pieces from a subcomm are not
27597929ea7SJunchao Zhang        necessarily consecutive in yg. That depends on how PETSC_COMM_WORLD is split. In our case, subcomm0 is made of rank
27697929ea7SJunchao Zhang        0, 3, 6 etc from PETSC_COMM_WORLD. So subcomm0's pieces are interleaved with pieces from other subcomms in yg.
27797929ea7SJunchao Zhang     */
27897929ea7SJunchao Zhang     ierr = VecGetArray(y,&yvalue);CHKERRQ(ierr);
27997929ea7SJunchao Zhang     ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,yvalue,&yg);CHKERRQ(ierr);
28097929ea7SJunchao Zhang     ierr = PetscObjectSetName((PetscObject)yg,"yg_on_subcomms");CHKERRQ(ierr); /* Give a name to view yg clearly */
28197929ea7SJunchao Zhang 
28297929ea7SJunchao Zhang     /* The following two lines are key. From xstart, we know where to pull entries from x. Note that we get xstart from y,
28397929ea7SJunchao Zhang        since first entry of y on this rank is from x[xstart]. From ystart, we know where ot put entries to yg.
28497929ea7SJunchao Zhang      */
28597929ea7SJunchao Zhang     ierr = VecGetOwnershipRange(y,&xstart,NULL);CHKERRQ(ierr);
28697929ea7SJunchao Zhang     ierr = VecGetOwnershipRange(yg,&ystart,NULL);CHKERRQ(ierr);
28797929ea7SJunchao Zhang 
28897929ea7SJunchao Zhang     ierr = ISCreateStride(PETSC_COMM_SELF,n,xstart,1,&ix);CHKERRQ(ierr);
28997929ea7SJunchao Zhang     ierr = ISCreateStride(PETSC_COMM_SELF,n,ystart,1,&iy);CHKERRQ(ierr);
29097929ea7SJunchao Zhang     ierr = VecScatterCreate(x,ix,yg,iy,&vscat);CHKERRQ(ierr);
29197929ea7SJunchao Zhang     ierr = VecScatterBegin(vscat,x,yg,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
29297929ea7SJunchao Zhang     ierr = VecScatterEnd(vscat,x,yg,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
29397929ea7SJunchao Zhang 
29497929ea7SJunchao Zhang     /* View yg on PETSC_COMM_WORLD before destroying it. We shall see the interleaving effect in output. */
29597929ea7SJunchao Zhang     ierr = VecView(yg,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
29697929ea7SJunchao Zhang     ierr = VecDestroy(&yg);CHKERRQ(ierr);
29797929ea7SJunchao Zhang 
29897929ea7SJunchao Zhang     /* Restory yvalue so that processes in subcomm can use y from now on. */
29997929ea7SJunchao Zhang     ierr = VecRestoreArray(y,&yvalue);CHKERRQ(ierr);
30097929ea7SJunchao Zhang     ierr = VecScale(y,3.0);CHKERRQ(ierr);
30197929ea7SJunchao Zhang 
30297929ea7SJunchao Zhang     ierr = ISDestroy(&ix);CHKERRQ(ierr); /* One can also destroy ix, iy immediately after VecScatterCreate() */
30397929ea7SJunchao Zhang     ierr = ISDestroy(&iy);CHKERRQ(ierr);
30497929ea7SJunchao Zhang     ierr = VecDestroy(&x);CHKERRQ(ierr);
30597929ea7SJunchao Zhang     ierr = VecDestroy(&y);CHKERRQ(ierr);
30697929ea7SJunchao Zhang     ierr = VecScatterDestroy(&vscat);CHKERRQ(ierr);
30797929ea7SJunchao Zhang   } /* world2subs */
30897929ea7SJunchao Zhang 
309*ffc4695bSBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRMPI(ierr);
31097929ea7SJunchao Zhang   ierr = PetscFinalize();
31197929ea7SJunchao Zhang   return ierr;
31297929ea7SJunchao Zhang }
31397929ea7SJunchao Zhang 
31497929ea7SJunchao Zhang /*TEST
31597929ea7SJunchao Zhang 
31697929ea7SJunchao Zhang    build:
31797929ea7SJunchao Zhang      requires: !define(PETSC_HAVE_MPIUNI)
31897929ea7SJunchao Zhang    testset:
31997929ea7SJunchao Zhang      nsize: 7
32097929ea7SJunchao Zhang 
32197929ea7SJunchao Zhang      test:
32297929ea7SJunchao Zhang        suffix: 1
32397929ea7SJunchao Zhang        args: -world2sub
32497929ea7SJunchao Zhang 
32597929ea7SJunchao Zhang      test:
32697929ea7SJunchao Zhang        suffix: 2
32797929ea7SJunchao Zhang        args: -sub2sub
32897929ea7SJunchao Zhang 
32997929ea7SJunchao Zhang      test:
33097929ea7SJunchao Zhang        suffix: 3
33197929ea7SJunchao Zhang        args: -world2subs
33297929ea7SJunchao Zhang 
33397929ea7SJunchao Zhang      test:
33497929ea7SJunchao Zhang        suffix: 4
33597929ea7SJunchao Zhang        args: -world2sub -sf_type neighbor
33697929ea7SJunchao Zhang        output_file: output/ex9_1.out
33797929ea7SJunchao Zhang        # OpenMPI has a bug wrt MPI_Neighbor_alltoallv etc (https://github.com/open-mpi/ompi/pull/6782). Once the patch is in, we can remove !define(PETSC_HAVE_OMPI_MAJOR_VERSION)
33897929ea7SJunchao Zhang        requires: define(PETSC_HAVE_MPI_NEIGHBORHOOD_COLLECTIVES) !define(PETSC_HAVE_OMPI_MAJOR_VERSION)
33997929ea7SJunchao Zhang TEST*/
34097929ea7SJunchao Zhang 
341