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