197929ea7SJunchao Zhang 297929ea7SJunchao Zhang static char help[]= "Scatters from a parallel vector to a sequential vector. \n\ 397929ea7SJunchao Zhang uses block index sets\n\n"; 497929ea7SJunchao Zhang 597929ea7SJunchao Zhang #include <petscvec.h> 697929ea7SJunchao Zhang 797929ea7SJunchao Zhang int main(int argc,char **argv) 897929ea7SJunchao Zhang { 997929ea7SJunchao Zhang PetscInt bs=1,n=5,i,low; 1097929ea7SJunchao Zhang PetscInt ix0[3] = {5,7,9},iy0[3] = {1,2,4},ix1[3] = {2,3,4},iy1[3] = {0,1,3}; 1197929ea7SJunchao Zhang PetscMPIInt size,rank; 1297929ea7SJunchao Zhang PetscScalar *array; 1397929ea7SJunchao Zhang Vec x,y; 1497929ea7SJunchao Zhang IS isx,isy; 1597929ea7SJunchao Zhang VecScatter ctx; 1697929ea7SJunchao Zhang 17*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 185f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 195f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 2097929ea7SJunchao Zhang 212c71b3e2SJacob Faibussowitsch PetscCheckFalse(size <2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run more than one processor"); 2297929ea7SJunchao Zhang 235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL)); 2497929ea7SJunchao Zhang n = bs*n; 2597929ea7SJunchao Zhang 2697929ea7SJunchao Zhang /* Create vector x over shared memory */ 275f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(x,n,PETSC_DECIDE)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(x)); 3097929ea7SJunchao Zhang 315f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(x,&low,NULL)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(x,&array)); 3397929ea7SJunchao Zhang for (i=0; i<n; i++) { 3497929ea7SJunchao Zhang array[i] = (PetscScalar)(i + low); 3597929ea7SJunchao Zhang } 365f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(x,&array)); 3797929ea7SJunchao Zhang 3897929ea7SJunchao Zhang /* Create a sequential vector y */ 395f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&y)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(y,0.0)); 4197929ea7SJunchao Zhang 4297929ea7SJunchao Zhang /* Create two index sets */ 43dd400576SPatrick Sanan if (rank == 0) { 445f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy)); 4697929ea7SJunchao Zhang } else { 475f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy)); 4997929ea7SJunchao Zhang } 5097929ea7SJunchao Zhang 5197929ea7SJunchao Zhang if (rank == 10) { 525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(ISView(isx,PETSC_VIEWER_STDOUT_SELF)); 5497929ea7SJunchao Zhang } 5597929ea7SJunchao Zhang 565f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(x,isx,y,isy,&ctx)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterSetFromOptions(ctx)); 5897929ea7SJunchao Zhang 5997929ea7SJunchao Zhang /* Test forward vecscatter */ 605f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD)); 6297929ea7SJunchao Zhang if (rank == 0) { 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[%d] y:\n",rank)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_SELF)); 6597929ea7SJunchao Zhang } 6697929ea7SJunchao Zhang 6797929ea7SJunchao Zhang /* Test reverse vecscatter */ 685f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(y,-1.0)); 6997929ea7SJunchao Zhang if (rank) { 705f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(y,1.0/(size - 1))); 7197929ea7SJunchao Zhang } 7297929ea7SJunchao Zhang 735f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_REVERSE)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_REVERSE)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 7697929ea7SJunchao Zhang 7797929ea7SJunchao Zhang /* Free spaces */ 785f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&ctx)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isx)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isy)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 83*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 84*b122ec5aSJacob Faibussowitsch return 0; 8597929ea7SJunchao Zhang } 8697929ea7SJunchao Zhang 8797929ea7SJunchao Zhang /*TEST 8897929ea7SJunchao Zhang 8997929ea7SJunchao Zhang test: 9097929ea7SJunchao Zhang nsize: 3 9197929ea7SJunchao Zhang 9297929ea7SJunchao Zhang TEST*/ 93