1 2 static char help[]= "Scatters from a parallel vector to a sequential vector. \n\ 3 uses block index sets\n\n"; 4 5 #include <petscvec.h> 6 7 int main(int argc,char **argv) 8 { 9 PetscInt bs=1,n=5,i,low; 10 PetscInt ix0[3] = {5,7,9},iy0[3] = {1,2,4},ix1[3] = {2,3,4},iy1[3] = {0,1,3}; 11 PetscMPIInt size,rank; 12 PetscScalar *array; 13 Vec x,y; 14 IS isx,isy; 15 VecScatter ctx; 16 17 PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 18 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 19 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 20 21 PetscCheck(size >=2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run more than one processor"); 22 23 PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL)); 24 n = bs*n; 25 26 /* Create vector x over shared memory */ 27 PetscCall(VecCreate(PETSC_COMM_WORLD,&x)); 28 PetscCall(VecSetSizes(x,n,PETSC_DECIDE)); 29 PetscCall(VecSetFromOptions(x)); 30 31 PetscCall(VecGetOwnershipRange(x,&low,NULL)); 32 PetscCall(VecGetArray(x,&array)); 33 for (i=0; i<n; i++) { 34 array[i] = (PetscScalar)(i + low); 35 } 36 PetscCall(VecRestoreArray(x,&array)); 37 38 /* Create a sequential vector y */ 39 PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&y)); 40 PetscCall(VecSet(y,0.0)); 41 42 /* Create two index sets */ 43 if (rank == 0) { 44 PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx)); 45 PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy)); 46 } else { 47 PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx)); 48 PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy)); 49 } 50 51 if (rank == 10) { 52 PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank)); 53 PetscCall(ISView(isx,PETSC_VIEWER_STDOUT_SELF)); 54 } 55 56 PetscCall(VecScatterCreate(x,isx,y,isy,&ctx)); 57 PetscCall(VecScatterSetFromOptions(ctx)); 58 59 /* Test forward vecscatter */ 60 PetscCall(VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD)); 61 PetscCall(VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD)); 62 if (rank == 0) { 63 PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] y:\n",rank)); 64 PetscCall(VecView(y,PETSC_VIEWER_STDOUT_SELF)); 65 } 66 67 /* Test reverse vecscatter */ 68 PetscCall(VecScale(y,-1.0)); 69 if (rank) { 70 PetscCall(VecScale(y,1.0/(size - 1))); 71 } 72 73 PetscCall(VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_REVERSE)); 74 PetscCall(VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_REVERSE)); 75 PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 76 77 /* Free spaces */ 78 PetscCall(VecScatterDestroy(&ctx)); 79 PetscCall(ISDestroy(&isx)); 80 PetscCall(ISDestroy(&isy)); 81 PetscCall(VecDestroy(&x)); 82 PetscCall(VecDestroy(&y)); 83 PetscCall(PetscFinalize()); 84 return 0; 85 } 86 87 /*TEST 88 89 test: 90 nsize: 3 91 92 TEST*/ 93