1 static char help[]= " Test VecScatter with x, y on different communicators\n\n"; 2 3 #include <petscvec.h> 4 5 int main(int argc,char **argv) 6 { 7 PetscErrorCode ierr; 8 PetscInt i,n=5,rstart; 9 PetscScalar *val; 10 const PetscScalar *dat; 11 Vec x,y1,y2; 12 MPI_Comm newcomm; 13 PetscMPIInt nproc,rank; 14 IS ix; 15 VecScatter vscat1,vscat2; 16 17 PetscFunctionBegin; 18 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 19 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&nproc);CHKERRMPI(ierr); 20 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 21 22 if (nproc != 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This test must run with exactly two MPI ranks\n"); 23 24 /* Create MPI vectors x and y, which are on the same comm (i.e., MPI_IDENT) */ 25 ierr = VecCreateMPI(PETSC_COMM_WORLD,n,PETSC_DECIDE,&x);CHKERRQ(ierr); 26 ierr = VecDuplicate(x,&y1);CHKERRQ(ierr); 27 ierr = VecGetOwnershipRange(x,&rstart,NULL);CHKERRQ(ierr); 28 29 /* Set x's value locally. x would be {0., 1., 2., ..., 9.} */ 30 ierr = VecGetArray(x,&val);CHKERRQ(ierr); 31 for (i=0; i<n; i++) val[i] = rstart + i; 32 ierr = VecRestoreArray(x,&val);CHKERRQ(ierr); 33 34 /* Create index set ix = {0, 1, 2, ..., 9} */ 35 ierr = ISCreateStride(PETSC_COMM_WORLD,n,rstart,1,&ix);CHKERRQ(ierr); 36 37 /* Create newcomm that reverses processes in x's comm, and then create y2 on it*/ 38 ierr = MPI_Comm_split(PETSC_COMM_WORLD,0/*color*/,-rank/*key*/,&newcomm);CHKERRMPI(ierr); 39 ierr = VecCreateMPI(newcomm,n,PETSC_DECIDE,&y2);CHKERRQ(ierr); 40 41 /* It looks vscat1/2 are the same, but actually not. y2 is on a different communicator than x */ 42 ierr = VecScatterCreate(x,ix,y1,ix,&vscat1);CHKERRQ(ierr); 43 ierr = VecScatterCreate(x,ix,y2,ix,&vscat2);CHKERRQ(ierr); 44 45 ierr = VecScatterBegin(vscat1,x,y1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 46 ierr = VecScatterBegin(vscat2,x,y2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 47 ierr = VecScatterEnd (vscat1,x,y1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 48 ierr = VecScatterEnd (vscat2,x,y2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 49 50 /* View on rank 0 of x's comm, which is PETSC_COMM_WORLD */ 51 if (!rank) { 52 /* Print the part of x on rank 0, which is 0 1 2 3 4 */ 53 ierr = PetscPrintf(PETSC_COMM_SELF,"\nOn rank 0 of PETSC_COMM_WORLD, x = ");CHKERRQ(ierr); 54 ierr = VecGetArrayRead(x,&dat);CHKERRQ(ierr); 55 for (i=0; i<n; i++) {ierr = PetscPrintf(PETSC_COMM_SELF," %.0f",PetscRealPart(dat[i]));CHKERRQ(ierr);} 56 ierr = VecRestoreArrayRead(x,&dat);CHKERRQ(ierr); 57 58 /* Print the part of y1 on rank 0, which is 0 1 2 3 4 */ 59 ierr = PetscPrintf(PETSC_COMM_SELF,"\nOn rank 0 of PETSC_COMM_WORLD, y1 = ");CHKERRQ(ierr); 60 ierr = VecGetArrayRead(y1,&dat);CHKERRQ(ierr); 61 for (i=0; i<n; i++) {ierr = PetscPrintf(PETSC_COMM_SELF," %.0f",PetscRealPart(dat[i]));CHKERRQ(ierr);} 62 ierr = VecRestoreArrayRead(y1,&dat);CHKERRQ(ierr); 63 64 /* Print the part of y2 on rank 0, which is 5 6 7 8 9 since y2 swapped the processes of PETSC_COMM_WORLD */ 65 ierr = PetscPrintf(PETSC_COMM_SELF,"\nOn rank 0 of PETSC_COMM_WORLD, y2 = ");CHKERRQ(ierr); 66 ierr = VecGetArrayRead(y2,&dat);CHKERRQ(ierr); 67 for (i=0; i<n; i++) {ierr = PetscPrintf(PETSC_COMM_SELF," %.0f",PetscRealPart(dat[i]));CHKERRQ(ierr);} 68 ierr = VecRestoreArrayRead(y2,&dat);CHKERRQ(ierr); 69 ierr = PetscPrintf(PETSC_COMM_SELF,"\n");CHKERRQ(ierr); 70 } 71 72 ierr = ISDestroy(&ix);CHKERRQ(ierr); 73 ierr = VecDestroy(&x);CHKERRQ(ierr); 74 ierr = VecDestroy(&y1);CHKERRQ(ierr); 75 ierr = VecDestroy(&y2);CHKERRQ(ierr); 76 ierr = VecScatterDestroy(&vscat1);CHKERRQ(ierr); 77 ierr = VecScatterDestroy(&vscat2);CHKERRQ(ierr); 78 ierr = MPI_Comm_free(&newcomm);CHKERRMPI(ierr); 79 ierr = PetscFinalize(); 80 return ierr; 81 } 82 83 /*TEST 84 85 test: 86 nsize: 2 87 requires: double 88 TEST*/ 89 90