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 PetscInt i,bs,n,low,high; 797929ea7SJunchao Zhang PetscMPIInt nproc,rank; 897929ea7SJunchao Zhang Vec x,y,z; 997929ea7SJunchao Zhang IS ix,iy; 1097929ea7SJunchao Zhang VecScatter vscat; 1197929ea7SJunchao Zhang const PetscScalar *yv; 1297929ea7SJunchao Zhang 13*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 145f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&nproc)); 155f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 162c71b3e2SJacob Faibussowitsch PetscCheckFalse(nproc != 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This test can only run on two MPI ranks"); 1797929ea7SJunchao Zhang 1897929ea7SJunchao Zhang /* Create an MPI vector x of size 12 on two processes, and set x = {0, 1, 2, .., 11} */ 195f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,6,PETSC_DECIDE,&x)); 205f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(x,&low,&high)); 215f80ce2aSJacob Faibussowitsch for (i=low; i<high; i++) CHKERRQ(VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(x)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(x)); 2497929ea7SJunchao Zhang 2597929ea7SJunchao Zhang /* Create a seq vector y, and a parallel to sequential (PtoS) vecscatter to scatter x to y */ 26dd400576SPatrick Sanan if (rank == 0) { 2797929ea7SJunchao 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 */ 2897929ea7SJunchao Zhang PetscInt idx[2]={0,2}; 2997929ea7SJunchao Zhang PetscInt idy[2]={0,1}; 3097929ea7SJunchao Zhang n = 6; 3197929ea7SJunchao Zhang bs = 3; 325f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&y)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,2,idx,PETSC_COPY_VALUES,&ix)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,2,idy,PETSC_COPY_VALUES,&iy)); 3597929ea7SJunchao Zhang } else { 3697929ea7SJunchao 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 */ 3797929ea7SJunchao Zhang PetscInt idx[2]= {2,5}; 3897929ea7SJunchao Zhang PetscInt idy[2]= {0,1}; 3997929ea7SJunchao Zhang n = 4; 4097929ea7SJunchao Zhang bs = 2; 415f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&y)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,2,idx,PETSC_COPY_VALUES,&ix)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,2,idy,PETSC_COPY_VALUES,&iy)); 4497929ea7SJunchao Zhang } 455f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(x,ix,y,iy,&vscat)); 4697929ea7SJunchao Zhang 4797929ea7SJunchao Zhang /* Do the vecscatter */ 485f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD)); 5097929ea7SJunchao Zhang 5197929ea7SJunchao Zhang /* Print y. Since y is sequential, we put y in a parallel z to print its value on both ranks */ 525f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(y,&yv)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,yv,&z)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(z,PETSC_VIEWER_STDOUT_WORLD)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(y,&yv)); 5697929ea7SJunchao Zhang 575f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&ix)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iy)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&z)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&vscat)); 6397929ea7SJunchao Zhang 64*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 65*b122ec5aSJacob Faibussowitsch return 0; 6697929ea7SJunchao Zhang } 6797929ea7SJunchao Zhang 6897929ea7SJunchao Zhang /*TEST 6997929ea7SJunchao Zhang 7097929ea7SJunchao Zhang test: 7197929ea7SJunchao Zhang nsize: 2 7297929ea7SJunchao Zhang args: 7397929ea7SJunchao Zhang requires: 7497929ea7SJunchao Zhang TEST*/ 75