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