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