197929ea7SJunchao Zhang 297929ea7SJunchao Zhang static char help[]= "Scatters from a sequential vector to a parallel 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 PetscViewer sviewer; 1797929ea7SJunchao Zhang 189566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 2197929ea7SJunchao Zhang 22*08401ef6SPierre Jolivet PetscCheck(size >=2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run more than one processor"); 2397929ea7SJunchao Zhang 249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL)); 2597929ea7SJunchao Zhang n = bs*n; 2697929ea7SJunchao Zhang 2797929ea7SJunchao Zhang /* Create vector x over shared memory */ 289566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&x)); 299566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x,n,PETSC_DECIDE)); 309566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 3197929ea7SJunchao Zhang 329566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x,&low,NULL)); 339566063dSJacob Faibussowitsch PetscCall(VecGetArray(x,&array)); 3497929ea7SJunchao Zhang for (i=0; i<n; i++) { 3597929ea7SJunchao Zhang array[i] = (PetscScalar)(i + low); 3697929ea7SJunchao Zhang } 379566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x,&array)); 3897929ea7SJunchao Zhang 3997929ea7SJunchao Zhang /* Create a sequential vector y */ 409566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&y)); 419566063dSJacob Faibussowitsch PetscCall(VecGetArray(y,&array)); 4297929ea7SJunchao Zhang for (i=0; i<n; i++) { 4397929ea7SJunchao Zhang array[i] = (PetscScalar)(i + 100*rank); 4497929ea7SJunchao Zhang } 459566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y,&array)); 4697929ea7SJunchao Zhang 4797929ea7SJunchao Zhang /* Create two index sets */ 48dd400576SPatrick Sanan if (rank == 0) { 499566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx)); 509566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy)); 5197929ea7SJunchao Zhang } else { 529566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx)); 539566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy)); 5497929ea7SJunchao Zhang } 5597929ea7SJunchao Zhang 5697929ea7SJunchao Zhang if (rank == 10) { 579566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank)); 589566063dSJacob Faibussowitsch PetscCall(ISView(isx,PETSC_VIEWER_STDOUT_SELF)); 5997929ea7SJunchao Zhang } 6097929ea7SJunchao Zhang 619566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(y,isy,x,isx,&ctx)); 629566063dSJacob Faibussowitsch PetscCall(VecScatterSetFromOptions(ctx)); 6397929ea7SJunchao Zhang 6497929ea7SJunchao Zhang /* Test forward vecscatter */ 659566063dSJacob Faibussowitsch PetscCall(VecSet(x,0.0)); 669566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_FORWARD)); 679566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_FORWARD)); 689566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 6997929ea7SJunchao Zhang 7097929ea7SJunchao Zhang /* Test reverse vecscatter */ 719566063dSJacob Faibussowitsch PetscCall(VecScale(x,-1.0)); 729566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_REVERSE)); 739566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_REVERSE)); 749566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 7597929ea7SJunchao Zhang if (rank == 1) { 769566063dSJacob Faibussowitsch PetscCall(VecView(y,sviewer)); 7797929ea7SJunchao Zhang } 789566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 7997929ea7SJunchao Zhang 8097929ea7SJunchao Zhang /* Free spaces */ 819566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ctx)); 829566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isx)); 839566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isy)); 849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 869566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 87b122ec5aSJacob Faibussowitsch return 0; 8897929ea7SJunchao Zhang } 8997929ea7SJunchao Zhang 9097929ea7SJunchao Zhang /*TEST 9197929ea7SJunchao Zhang 9297929ea7SJunchao Zhang test: 9397929ea7SJunchao Zhang nsize: 3 9497929ea7SJunchao Zhang 9597929ea7SJunchao Zhang TEST*/ 96