197929ea7SJunchao Zhang 297929ea7SJunchao Zhang static char help[]= "Scatters between parallel vectors. \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,N,i,low; 1097929ea7SJunchao Zhang PetscInt ix0[3] = {5,7,9},iy0[3] = {1,2,4},ix1[3] = {2,3,1},iy1[3] = {0,3,9}; 1197929ea7SJunchao Zhang PetscMPIInt size,rank; 1297929ea7SJunchao Zhang PetscScalar *array; 1397929ea7SJunchao Zhang Vec x,x1,y; 1497929ea7SJunchao Zhang IS isx,isy; 1597929ea7SJunchao Zhang VecScatter ctx; 1697929ea7SJunchao Zhang VecScatterType type; 1797929ea7SJunchao Zhang PetscBool flg; 1897929ea7SJunchao Zhang 19*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 20*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 21*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 2297929ea7SJunchao Zhang 232c71b3e2SJacob Faibussowitsch PetscCheckFalse(size <2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run more than one processor"); 2497929ea7SJunchao Zhang 25*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL)); 2697929ea7SJunchao Zhang n = bs*n; 2797929ea7SJunchao Zhang 2897929ea7SJunchao Zhang /* Create vector x over shared memory */ 29*9566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&x)); 30*9566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x,n,PETSC_DECIDE)); 31*9566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 3297929ea7SJunchao Zhang 33*9566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x,&low,NULL)); 34*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(x,&array)); 3597929ea7SJunchao Zhang for (i=0; i<n; i++) { 3697929ea7SJunchao Zhang array[i] = (PetscScalar)(i + low); 3797929ea7SJunchao Zhang } 38*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x,&array)); 39*9566063dSJacob Faibussowitsch /* PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); */ 4097929ea7SJunchao Zhang 4197929ea7SJunchao Zhang /* Test some vector functions */ 42*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x)); 43*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x)); 4497929ea7SJunchao Zhang 45*9566063dSJacob Faibussowitsch PetscCall(VecGetSize(x,&N)); 46*9566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x,&n)); 4797929ea7SJunchao Zhang 48*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&x1)); 49*9566063dSJacob Faibussowitsch PetscCall(VecCopy(x,x1)); 50*9566063dSJacob Faibussowitsch PetscCall(VecEqual(x,x1,&flg)); 5128b400f6SJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_WRONG,"x1 != x"); 5297929ea7SJunchao Zhang 53*9566063dSJacob Faibussowitsch PetscCall(VecScale(x1,2.0)); 54*9566063dSJacob Faibussowitsch PetscCall(VecSet(x1,10.0)); 55*9566063dSJacob Faibussowitsch /* PetscCall(VecView(x1,PETSC_VIEWER_STDOUT_WORLD)); */ 5697929ea7SJunchao Zhang 5797929ea7SJunchao Zhang /* Create vector y over shared memory */ 58*9566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&y)); 59*9566063dSJacob Faibussowitsch PetscCall(VecSetSizes(y,n,PETSC_DECIDE)); 60*9566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(y)); 61*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(y,&array)); 6297929ea7SJunchao Zhang for (i=0; i<n; i++) { 6397929ea7SJunchao Zhang array[i] = -(PetscScalar) (i + 100*rank); 6497929ea7SJunchao Zhang } 65*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y,&array)); 66*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(y)); 67*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(y)); 68*9566063dSJacob Faibussowitsch /* PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); */ 6997929ea7SJunchao Zhang 7097929ea7SJunchao Zhang /* Create two index sets */ 71dd400576SPatrick Sanan if (rank == 0) { 72*9566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx)); 73*9566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy)); 7497929ea7SJunchao Zhang } else { 75*9566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx)); 76*9566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy)); 7797929ea7SJunchao Zhang } 7897929ea7SJunchao Zhang 7997929ea7SJunchao Zhang if (rank == 10) { 80*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank)); 81*9566063dSJacob Faibussowitsch PetscCall(ISView(isx,PETSC_VIEWER_STDOUT_SELF)); 82*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n[%d] isy:\n",rank)); 83*9566063dSJacob Faibussowitsch PetscCall(ISView(isy,PETSC_VIEWER_STDOUT_SELF)); 8497929ea7SJunchao Zhang } 8597929ea7SJunchao Zhang 8697929ea7SJunchao Zhang /* Create Vector scatter */ 87*9566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x,isx,y,isy,&ctx)); 88*9566063dSJacob Faibussowitsch PetscCall(VecScatterSetFromOptions(ctx)); 89*9566063dSJacob Faibussowitsch PetscCall(VecScatterGetType(ctx,&type)); 90*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"scatter type %s\n",type)); 9197929ea7SJunchao Zhang 9297929ea7SJunchao Zhang /* Test forward vecscatter */ 93*9566063dSJacob Faibussowitsch PetscCall(VecSet(y,0.0)); 94*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD)); 95*9566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD)); 96*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSCATTER_FORWARD y:\n")); 97*9566063dSJacob Faibussowitsch PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); 9897929ea7SJunchao Zhang 9997929ea7SJunchao Zhang /* Test reverse vecscatter */ 100*9566063dSJacob Faibussowitsch PetscCall(VecSet(x,0.0)); 101*9566063dSJacob Faibussowitsch PetscCall(VecSet(y,0.0)); 102*9566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(y,&low,NULL)); 103*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(y,&array)); 10497929ea7SJunchao Zhang for (i=0; i<n; i++) { 10597929ea7SJunchao Zhang array[i] = (PetscScalar)(i+ low); 10697929ea7SJunchao Zhang } 107*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y,&array)); 108*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_REVERSE)); 109*9566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_REVERSE)); 110*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSCATTER_REVERSE x:\n")); 111*9566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 11297929ea7SJunchao Zhang 11397929ea7SJunchao Zhang /* Free objects */ 114*9566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ctx)); 115*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isx)); 116*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isy)); 117*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 118*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x1)); 119*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 120*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 121b122ec5aSJacob Faibussowitsch return 0; 12297929ea7SJunchao Zhang } 12397929ea7SJunchao Zhang 12497929ea7SJunchao Zhang /*TEST 12597929ea7SJunchao Zhang 12697929ea7SJunchao Zhang test: 12797929ea7SJunchao Zhang nsize: 2 12897929ea7SJunchao Zhang 12997929ea7SJunchao Zhang TEST*/ 130