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