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 PetscErrorCode ierr; 897929ea7SJunchao Zhang PetscInt i,n=5,rstart; 997929ea7SJunchao Zhang PetscScalar *val; 1097929ea7SJunchao Zhang const PetscScalar *dat; 1197929ea7SJunchao Zhang Vec x,y1,y2; 1297929ea7SJunchao Zhang MPI_Comm newcomm; 1397929ea7SJunchao Zhang PetscMPIInt nproc,rank; 1497929ea7SJunchao Zhang IS ix; 1597929ea7SJunchao Zhang VecScatter vscat1,vscat2; 1697929ea7SJunchao Zhang 1797929ea7SJunchao Zhang PetscFunctionBegin; 1897929ea7SJunchao Zhang ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 19*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&nproc)); 20*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 2197929ea7SJunchao Zhang 222c71b3e2SJacob Faibussowitsch PetscCheckFalse(nproc != 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This test must run with exactly two MPI ranks"); 2397929ea7SJunchao Zhang 2497929ea7SJunchao Zhang /* Create MPI vectors x and y, which are on the same comm (i.e., MPI_IDENT) */ 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,n,PETSC_DECIDE,&x)); 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&y1)); 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(x,&rstart,NULL)); 2897929ea7SJunchao Zhang 2997929ea7SJunchao Zhang /* Set x's value locally. x would be {0., 1., 2., ..., 9.} */ 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(x,&val)); 3197929ea7SJunchao Zhang for (i=0; i<n; i++) val[i] = rstart + i; 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(x,&val)); 3397929ea7SJunchao Zhang 3497929ea7SJunchao Zhang /* Create index set ix = {0, 1, 2, ..., 9} */ 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_WORLD,n,rstart,1,&ix)); 3697929ea7SJunchao Zhang 3797929ea7SJunchao Zhang /* Create newcomm that reverses processes in x's comm, and then create y2 on it*/ 38*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_split(PETSC_COMM_WORLD,0/*color*/,-rank/*key*/,&newcomm)); 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(newcomm,n,PETSC_DECIDE,&y2)); 4097929ea7SJunchao Zhang 4197929ea7SJunchao Zhang /* It looks vscat1/2 are the same, but actually not. y2 is on a different communicator than x */ 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(x,ix,y1,ix,&vscat1)); 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(x,ix,y2,ix,&vscat2)); 4497929ea7SJunchao Zhang 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(vscat1,x,y1,INSERT_VALUES,SCATTER_FORWARD)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(vscat2,x,y2,INSERT_VALUES,SCATTER_FORWARD)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd (vscat1,x,y1,INSERT_VALUES,SCATTER_FORWARD)); 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd (vscat2,x,y2,INSERT_VALUES,SCATTER_FORWARD)); 4997929ea7SJunchao Zhang 5097929ea7SJunchao Zhang /* View on rank 0 of x's comm, which is PETSC_COMM_WORLD */ 51dd400576SPatrick Sanan if (rank == 0) { 5297929ea7SJunchao Zhang /* Print the part of x on rank 0, which is 0 1 2 3 4 */ 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\nOn rank 0 of PETSC_COMM_WORLD, x = ")); 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&dat)); 55*5f80ce2aSJacob Faibussowitsch for (i=0; i<n; i++) CHKERRQ(PetscPrintf(PETSC_COMM_SELF," %.0g",(double)PetscRealPart(dat[i]))); 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&dat)); 5797929ea7SJunchao Zhang 5897929ea7SJunchao Zhang /* Print the part of y1 on rank 0, which is 0 1 2 3 4 */ 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\nOn rank 0 of PETSC_COMM_WORLD, y1 = ")); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(y1,&dat)); 61*5f80ce2aSJacob Faibussowitsch for (i=0; i<n; i++) CHKERRQ(PetscPrintf(PETSC_COMM_SELF," %.0g",(double)PetscRealPart(dat[i]))); 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(y1,&dat)); 6397929ea7SJunchao Zhang 6497929ea7SJunchao 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 */ 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\nOn rank 0 of PETSC_COMM_WORLD, y2 = ")); 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(y2,&dat)); 67*5f80ce2aSJacob Faibussowitsch for (i=0; i<n; i++) CHKERRQ(PetscPrintf(PETSC_COMM_SELF," %.0g",(double)PetscRealPart(dat[i]))); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(y2,&dat)); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n")); 7097929ea7SJunchao Zhang } 7197929ea7SJunchao Zhang 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&ix)); 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y1)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y2)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&vscat1)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&vscat2)); 78*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_free(&newcomm)); 7997929ea7SJunchao Zhang ierr = PetscFinalize(); 8097929ea7SJunchao Zhang return ierr; 8197929ea7SJunchao Zhang } 8297929ea7SJunchao Zhang 8397929ea7SJunchao Zhang /*TEST 8497929ea7SJunchao Zhang 8597929ea7SJunchao Zhang test: 8697929ea7SJunchao Zhang nsize: 2 8797929ea7SJunchao Zhang requires: double 8897929ea7SJunchao Zhang TEST*/ 89