xref: /petsc/src/vec/is/sf/tests/ex7.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
197929ea7SJunchao Zhang static char help[]= "Test vecscatter of different block sizes across processes\n\n";
297929ea7SJunchao Zhang 
397929ea7SJunchao Zhang #include <petscvec.h>
497929ea7SJunchao Zhang int main(int argc,char **argv)
597929ea7SJunchao Zhang {
697929ea7SJunchao Zhang   PetscErrorCode     ierr;
797929ea7SJunchao Zhang   PetscInt           i,bs,n,low,high;
897929ea7SJunchao Zhang   PetscMPIInt        nproc,rank;
997929ea7SJunchao Zhang   Vec                x,y,z;
1097929ea7SJunchao Zhang   IS                 ix,iy;
1197929ea7SJunchao Zhang   VecScatter         vscat;
1297929ea7SJunchao Zhang   const PetscScalar  *yv;
1397929ea7SJunchao Zhang 
1497929ea7SJunchao Zhang   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
15*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&nproc));
16*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
172c71b3e2SJacob Faibussowitsch   PetscCheckFalse(nproc != 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This test can only run on two MPI ranks");
1897929ea7SJunchao Zhang 
1997929ea7SJunchao Zhang   /* Create an MPI vector x of size 12 on two processes, and set x = {0, 1, 2, .., 11} */
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,6,PETSC_DECIDE,&x));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(x,&low,&high));
22*5f80ce2aSJacob Faibussowitsch   for (i=low; i<high; i++) CHKERRQ(VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES));
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(x));
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(x));
2597929ea7SJunchao Zhang 
2697929ea7SJunchao Zhang   /* Create a seq vector y, and a parallel to sequential (PtoS) vecscatter to scatter x to y */
27dd400576SPatrick Sanan   if (rank == 0) {
2897929ea7SJunchao Zhang     /* On rank 0, seq y is of size 6. We will scatter x[0,1,2,6,7,8] to y[0,1,2,3,4,5] using IS with bs=3 */
2997929ea7SJunchao Zhang     PetscInt idx[2]={0,2};
3097929ea7SJunchao Zhang     PetscInt idy[2]={0,1};
3197929ea7SJunchao Zhang     n    = 6;
3297929ea7SJunchao Zhang     bs   = 3;
33*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&y));
34*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,2,idx,PETSC_COPY_VALUES,&ix));
35*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,2,idy,PETSC_COPY_VALUES,&iy));
3697929ea7SJunchao Zhang   } else {
3797929ea7SJunchao Zhang     /* On rank 1, seq y is of size 4. We will scatter x[4,5,10,11] to y[0,1,2,3] using IS with bs=2 */
3897929ea7SJunchao Zhang     PetscInt idx[2]= {2,5};
3997929ea7SJunchao Zhang     PetscInt idy[2]= {0,1};
4097929ea7SJunchao Zhang     n    = 4;
4197929ea7SJunchao Zhang     bs   = 2;
42*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&y));
43*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,2,idx,PETSC_COPY_VALUES,&ix));
44*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,2,idy,PETSC_COPY_VALUES,&iy));
4597929ea7SJunchao Zhang   }
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(x,ix,y,iy,&vscat));
4797929ea7SJunchao Zhang 
4897929ea7SJunchao Zhang   /* Do the vecscatter */
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
5197929ea7SJunchao Zhang 
5297929ea7SJunchao Zhang   /* Print y. Since y is sequential, we put y in a parallel z to print its value on both ranks */
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(y,&yv));
54*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,yv,&z));
55*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(z,PETSC_VIEWER_STDOUT_WORLD));
56*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(y,&yv));
5797929ea7SJunchao Zhang 
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&ix));
59*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&iy));
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&z));
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&vscat));
6497929ea7SJunchao Zhang 
6597929ea7SJunchao Zhang   ierr = PetscFinalize();
6697929ea7SJunchao Zhang   return ierr;
6797929ea7SJunchao Zhang }
6897929ea7SJunchao Zhang 
6997929ea7SJunchao Zhang /*TEST
7097929ea7SJunchao Zhang 
7197929ea7SJunchao Zhang    test:
7297929ea7SJunchao Zhang       nsize: 2
7397929ea7SJunchao Zhang       args:
7497929ea7SJunchao Zhang       requires:
7597929ea7SJunchao Zhang TEST*/
76