197929ea7SJunchao Zhang static char help[]= " Test VecScatter with x, y on different communicators\n\n"; 297929ea7SJunchao Zhang 397929ea7SJunchao Zhang #include <petscvec.h> 497929ea7SJunchao Zhang 597929ea7SJunchao Zhang int main(int argc,char **argv) 697929ea7SJunchao Zhang { 797929ea7SJunchao Zhang PetscInt i,n=5,rstart; 897929ea7SJunchao Zhang PetscScalar *val; 997929ea7SJunchao Zhang const PetscScalar *dat; 1097929ea7SJunchao Zhang Vec x,y1,y2; 1197929ea7SJunchao Zhang MPI_Comm newcomm; 1297929ea7SJunchao Zhang PetscMPIInt nproc,rank; 1397929ea7SJunchao Zhang IS ix; 1497929ea7SJunchao Zhang VecScatter vscat1,vscat2; 1597929ea7SJunchao Zhang 1697929ea7SJunchao Zhang PetscFunctionBegin; 179566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&nproc)); 199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 2097929ea7SJunchao Zhang 21*08401ef6SPierre Jolivet PetscCheck(nproc == 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This test must run with exactly two MPI ranks"); 2297929ea7SJunchao Zhang 2397929ea7SJunchao Zhang /* Create MPI vectors x and y, which are on the same comm (i.e., MPI_IDENT) */ 249566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PETSC_COMM_WORLD,n,PETSC_DECIDE,&x)); 259566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&y1)); 269566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x,&rstart,NULL)); 2797929ea7SJunchao Zhang 2897929ea7SJunchao Zhang /* Set x's value locally. x would be {0., 1., 2., ..., 9.} */ 299566063dSJacob Faibussowitsch PetscCall(VecGetArray(x,&val)); 3097929ea7SJunchao Zhang for (i=0; i<n; i++) val[i] = rstart + i; 319566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x,&val)); 3297929ea7SJunchao Zhang 3397929ea7SJunchao Zhang /* Create index set ix = {0, 1, 2, ..., 9} */ 349566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,n,rstart,1,&ix)); 3597929ea7SJunchao Zhang 3697929ea7SJunchao Zhang /* Create newcomm that reverses processes in x's comm, and then create y2 on it*/ 379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_split(PETSC_COMM_WORLD,0/*color*/,-rank/*key*/,&newcomm)); 389566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(newcomm,n,PETSC_DECIDE,&y2)); 3997929ea7SJunchao Zhang 4097929ea7SJunchao Zhang /* It looks vscat1/2 are the same, but actually not. y2 is on a different communicator than x */ 419566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x,ix,y1,ix,&vscat1)); 429566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x,ix,y2,ix,&vscat2)); 4397929ea7SJunchao Zhang 449566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat1,x,y1,INSERT_VALUES,SCATTER_FORWARD)); 459566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat2,x,y2,INSERT_VALUES,SCATTER_FORWARD)); 469566063dSJacob Faibussowitsch PetscCall(VecScatterEnd (vscat1,x,y1,INSERT_VALUES,SCATTER_FORWARD)); 479566063dSJacob Faibussowitsch PetscCall(VecScatterEnd (vscat2,x,y2,INSERT_VALUES,SCATTER_FORWARD)); 4897929ea7SJunchao Zhang 4997929ea7SJunchao Zhang /* View on rank 0 of x's comm, which is PETSC_COMM_WORLD */ 50dd400576SPatrick Sanan if (rank == 0) { 5197929ea7SJunchao Zhang /* Print the part of x on rank 0, which is 0 1 2 3 4 */ 529566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\nOn rank 0 of PETSC_COMM_WORLD, x = ")); 539566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&dat)); 549566063dSJacob Faibussowitsch for (i=0; i<n; i++) PetscCall(PetscPrintf(PETSC_COMM_SELF," %.0g",(double)PetscRealPart(dat[i]))); 559566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&dat)); 5697929ea7SJunchao Zhang 5797929ea7SJunchao Zhang /* Print the part of y1 on rank 0, which is 0 1 2 3 4 */ 589566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\nOn rank 0 of PETSC_COMM_WORLD, y1 = ")); 599566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(y1,&dat)); 609566063dSJacob Faibussowitsch for (i=0; i<n; i++) PetscCall(PetscPrintf(PETSC_COMM_SELF," %.0g",(double)PetscRealPart(dat[i]))); 619566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(y1,&dat)); 6297929ea7SJunchao Zhang 6397929ea7SJunchao Zhang /* Print the part of y2 on rank 0, which is 5 6 7 8 9 since y2 swapped the processes of PETSC_COMM_WORLD */ 649566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\nOn rank 0 of PETSC_COMM_WORLD, y2 = ")); 659566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(y2,&dat)); 669566063dSJacob Faibussowitsch for (i=0; i<n; i++) PetscCall(PetscPrintf(PETSC_COMM_SELF," %.0g",(double)PetscRealPart(dat[i]))); 679566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(y2,&dat)); 689566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n")); 6997929ea7SJunchao Zhang } 7097929ea7SJunchao Zhang 719566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ix)); 729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y1)); 749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y2)); 759566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&vscat1)); 769566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&vscat2)); 779566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&newcomm)); 789566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 79b122ec5aSJacob Faibussowitsch return 0; 8097929ea7SJunchao Zhang } 8197929ea7SJunchao Zhang 8297929ea7SJunchao Zhang /*TEST 8397929ea7SJunchao Zhang 8497929ea7SJunchao Zhang test: 8597929ea7SJunchao Zhang nsize: 2 8697929ea7SJunchao Zhang requires: double 8797929ea7SJunchao Zhang TEST*/ 88