197929ea7SJunchao Zhang static char help[]= "Test vecscatter of different block sizes across processes\n\n"; 297929ea7SJunchao Zhang 397929ea7SJunchao Zhang #include <petscvec.h> 497929ea7SJunchao Zhang int main(int argc,char **argv) 597929ea7SJunchao Zhang { 697929ea7SJunchao Zhang PetscErrorCode ierr; 797929ea7SJunchao Zhang PetscInt i,bs,n,low,high; 897929ea7SJunchao Zhang PetscMPIInt nproc,rank; 997929ea7SJunchao Zhang Vec x,y,z; 1097929ea7SJunchao Zhang IS ix,iy; 1197929ea7SJunchao Zhang VecScatter vscat; 1297929ea7SJunchao Zhang const PetscScalar *yv; 1397929ea7SJunchao Zhang 1497929ea7SJunchao Zhang ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 15ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&nproc);CHKERRMPI(ierr); 16ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 17*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(nproc != 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This test can only run on two MPI ranks"); 1897929ea7SJunchao Zhang 1997929ea7SJunchao Zhang /* Create an MPI vector x of size 12 on two processes, and set x = {0, 1, 2, .., 11} */ 2097929ea7SJunchao Zhang ierr = VecCreateMPI(PETSC_COMM_WORLD,6,PETSC_DECIDE,&x);CHKERRQ(ierr); 2197929ea7SJunchao Zhang ierr = VecGetOwnershipRange(x,&low,&high);CHKERRQ(ierr); 2297929ea7SJunchao Zhang for (i=low; i<high; i++) {ierr = VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES);CHKERRQ(ierr);} 2397929ea7SJunchao Zhang ierr = VecAssemblyBegin(x);CHKERRQ(ierr); 2497929ea7SJunchao Zhang ierr = VecAssemblyEnd(x);CHKERRQ(ierr); 2597929ea7SJunchao Zhang 2697929ea7SJunchao Zhang /* Create a seq vector y, and a parallel to sequential (PtoS) vecscatter to scatter x to y */ 27dd400576SPatrick Sanan if (rank == 0) { 2897929ea7SJunchao Zhang /* On rank 0, seq y is of size 6. We will scatter x[0,1,2,6,7,8] to y[0,1,2,3,4,5] using IS with bs=3 */ 2997929ea7SJunchao Zhang PetscInt idx[2]={0,2}; 3097929ea7SJunchao Zhang PetscInt idy[2]={0,1}; 3197929ea7SJunchao Zhang n = 6; 3297929ea7SJunchao Zhang bs = 3; 3397929ea7SJunchao Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,n,&y);CHKERRQ(ierr); 3497929ea7SJunchao Zhang ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2,idx,PETSC_COPY_VALUES,&ix);CHKERRQ(ierr); 3597929ea7SJunchao Zhang ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2,idy,PETSC_COPY_VALUES,&iy);CHKERRQ(ierr); 3697929ea7SJunchao Zhang } else { 3797929ea7SJunchao Zhang /* On rank 1, seq y is of size 4. We will scatter x[4,5,10,11] to y[0,1,2,3] using IS with bs=2 */ 3897929ea7SJunchao Zhang PetscInt idx[2]= {2,5}; 3997929ea7SJunchao Zhang PetscInt idy[2]= {0,1}; 4097929ea7SJunchao Zhang n = 4; 4197929ea7SJunchao Zhang bs = 2; 4297929ea7SJunchao Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,n,&y);CHKERRQ(ierr); 4397929ea7SJunchao Zhang ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2,idx,PETSC_COPY_VALUES,&ix);CHKERRQ(ierr); 4497929ea7SJunchao Zhang ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2,idy,PETSC_COPY_VALUES,&iy);CHKERRQ(ierr); 4597929ea7SJunchao Zhang } 4697929ea7SJunchao Zhang ierr = VecScatterCreate(x,ix,y,iy,&vscat);CHKERRQ(ierr); 4797929ea7SJunchao Zhang 4897929ea7SJunchao Zhang /* Do the vecscatter */ 4997929ea7SJunchao Zhang ierr = VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5097929ea7SJunchao Zhang ierr = VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5197929ea7SJunchao Zhang 5297929ea7SJunchao Zhang /* Print y. Since y is sequential, we put y in a parallel z to print its value on both ranks */ 5397929ea7SJunchao Zhang ierr = VecGetArrayRead(y,&yv);CHKERRQ(ierr); 5497929ea7SJunchao Zhang ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,yv,&z);CHKERRQ(ierr); 5597929ea7SJunchao Zhang ierr = VecView(z,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 5697929ea7SJunchao Zhang ierr = VecRestoreArrayRead(y,&yv);CHKERRQ(ierr); 5797929ea7SJunchao Zhang 5897929ea7SJunchao Zhang ierr = ISDestroy(&ix);CHKERRQ(ierr); 5997929ea7SJunchao Zhang ierr = ISDestroy(&iy);CHKERRQ(ierr); 6097929ea7SJunchao Zhang ierr = VecDestroy(&x);CHKERRQ(ierr); 6197929ea7SJunchao Zhang ierr = VecDestroy(&y);CHKERRQ(ierr); 6297929ea7SJunchao Zhang ierr = VecDestroy(&z);CHKERRQ(ierr); 6397929ea7SJunchao Zhang ierr = VecScatterDestroy(&vscat);CHKERRQ(ierr); 6497929ea7SJunchao Zhang 6597929ea7SJunchao Zhang ierr = PetscFinalize(); 6697929ea7SJunchao Zhang return ierr; 6797929ea7SJunchao Zhang } 6897929ea7SJunchao Zhang 6997929ea7SJunchao Zhang /*TEST 7097929ea7SJunchao Zhang 7197929ea7SJunchao Zhang test: 7297929ea7SJunchao Zhang nsize: 2 7397929ea7SJunchao Zhang args: 7497929ea7SJunchao Zhang requires: 7597929ea7SJunchao Zhang TEST*/ 7697929ea7SJunchao Zhang 77