197929ea7SJunchao Zhang 297929ea7SJunchao Zhang static char help[]= "Scatters between parallel vectors. \n\ 397929ea7SJunchao Zhang uses block index sets\n\n"; 497929ea7SJunchao Zhang 597929ea7SJunchao Zhang #include <petscvec.h> 697929ea7SJunchao Zhang 797929ea7SJunchao Zhang int main(int argc,char **argv) 897929ea7SJunchao Zhang { 997929ea7SJunchao Zhang PetscErrorCode ierr; 1097929ea7SJunchao Zhang PetscInt bs=1,n=5,N,i,low; 1197929ea7SJunchao Zhang PetscInt ix0[3] = {5,7,9},iy0[3] = {1,2,4},ix1[3] = {2,3,1},iy1[3] = {0,3,9}; 1297929ea7SJunchao Zhang PetscMPIInt size,rank; 1397929ea7SJunchao Zhang PetscScalar *array; 1497929ea7SJunchao Zhang Vec x,x1,y; 1597929ea7SJunchao Zhang IS isx,isy; 1697929ea7SJunchao Zhang VecScatter ctx; 1797929ea7SJunchao Zhang VecScatterType type; 1897929ea7SJunchao Zhang PetscBool flg; 1997929ea7SJunchao Zhang 2097929ea7SJunchao Zhang ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 21ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 22ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 2397929ea7SJunchao Zhang 24*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(size <2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run more than one processor"); 2597929ea7SJunchao Zhang 2697929ea7SJunchao Zhang ierr = PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);CHKERRQ(ierr); 2797929ea7SJunchao Zhang n = bs*n; 2897929ea7SJunchao Zhang 2997929ea7SJunchao Zhang /* Create vector x over shared memory */ 3097929ea7SJunchao Zhang ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 3197929ea7SJunchao Zhang ierr = VecSetSizes(x,n,PETSC_DECIDE);CHKERRQ(ierr); 3297929ea7SJunchao Zhang ierr = VecSetFromOptions(x);CHKERRQ(ierr); 3397929ea7SJunchao Zhang 3497929ea7SJunchao Zhang ierr = VecGetOwnershipRange(x,&low,NULL);CHKERRQ(ierr); 3597929ea7SJunchao Zhang ierr = VecGetArray(x,&array);CHKERRQ(ierr); 3697929ea7SJunchao Zhang for (i=0; i<n; i++) { 3797929ea7SJunchao Zhang array[i] = (PetscScalar)(i + low); 3897929ea7SJunchao Zhang } 3997929ea7SJunchao Zhang ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 4097929ea7SJunchao Zhang /* ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 4197929ea7SJunchao Zhang 4297929ea7SJunchao Zhang /* Test some vector functions */ 4397929ea7SJunchao Zhang ierr = VecAssemblyBegin(x);CHKERRQ(ierr); 4497929ea7SJunchao Zhang ierr = VecAssemblyEnd(x);CHKERRQ(ierr); 4597929ea7SJunchao Zhang 4697929ea7SJunchao Zhang ierr = VecGetSize(x,&N);CHKERRQ(ierr); 4797929ea7SJunchao Zhang ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr); 4897929ea7SJunchao Zhang 4997929ea7SJunchao Zhang ierr = VecDuplicate(x,&x1);CHKERRQ(ierr); 5097929ea7SJunchao Zhang ierr = VecCopy(x,x1);CHKERRQ(ierr); 5197929ea7SJunchao Zhang ierr = VecEqual(x,x1,&flg);CHKERRQ(ierr); 52*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_WRONG,"x1 != x"); 5397929ea7SJunchao Zhang 5497929ea7SJunchao Zhang ierr = VecScale(x1,2.0);CHKERRQ(ierr); 5597929ea7SJunchao Zhang ierr = VecSet(x1,10.0);CHKERRQ(ierr); 5697929ea7SJunchao Zhang /* ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 5797929ea7SJunchao Zhang 5897929ea7SJunchao Zhang /* Create vector y over shared memory */ 5997929ea7SJunchao Zhang ierr = VecCreate(PETSC_COMM_WORLD,&y);CHKERRQ(ierr); 6097929ea7SJunchao Zhang ierr = VecSetSizes(y,n,PETSC_DECIDE);CHKERRQ(ierr); 6197929ea7SJunchao Zhang ierr = VecSetFromOptions(y);CHKERRQ(ierr); 6297929ea7SJunchao Zhang ierr = VecGetArray(y,&array);CHKERRQ(ierr); 6397929ea7SJunchao Zhang for (i=0; i<n; i++) { 6497929ea7SJunchao Zhang array[i] = -(PetscScalar) (i + 100*rank); 6597929ea7SJunchao Zhang } 6697929ea7SJunchao Zhang ierr = VecRestoreArray(y,&array);CHKERRQ(ierr); 6797929ea7SJunchao Zhang ierr = VecAssemblyBegin(y);CHKERRQ(ierr); 6897929ea7SJunchao Zhang ierr = VecAssemblyEnd(y);CHKERRQ(ierr); 6997929ea7SJunchao Zhang /* ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 7097929ea7SJunchao Zhang 7197929ea7SJunchao Zhang /* Create two index sets */ 72dd400576SPatrick Sanan if (rank == 0) { 7397929ea7SJunchao Zhang ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr); 7497929ea7SJunchao Zhang ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy);CHKERRQ(ierr); 7597929ea7SJunchao Zhang } else { 7697929ea7SJunchao Zhang ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr); 7797929ea7SJunchao Zhang ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy);CHKERRQ(ierr); 7897929ea7SJunchao Zhang } 7997929ea7SJunchao Zhang 8097929ea7SJunchao Zhang if (rank == 10) { 8197929ea7SJunchao Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank);CHKERRQ(ierr); 8297929ea7SJunchao Zhang ierr = ISView(isx,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 8397929ea7SJunchao Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"\n[%d] isy:\n",rank);CHKERRQ(ierr); 8497929ea7SJunchao Zhang ierr = ISView(isy,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 8597929ea7SJunchao Zhang } 8697929ea7SJunchao Zhang 8797929ea7SJunchao Zhang /* Create Vector scatter */ 8897929ea7SJunchao Zhang ierr = VecScatterCreate(x,isx,y,isy,&ctx);CHKERRQ(ierr); 8997929ea7SJunchao Zhang ierr = VecScatterSetFromOptions(ctx);CHKERRQ(ierr); 9097929ea7SJunchao Zhang ierr = VecScatterGetType(ctx,&type);CHKERRQ(ierr); 9197929ea7SJunchao Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"scatter type %s\n",type);CHKERRQ(ierr); 9297929ea7SJunchao Zhang 9397929ea7SJunchao Zhang /* Test forward vecscatter */ 9497929ea7SJunchao Zhang ierr = VecSet(y,0.0);CHKERRQ(ierr); 9597929ea7SJunchao Zhang ierr = VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9697929ea7SJunchao Zhang ierr = VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9797929ea7SJunchao Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\nSCATTER_FORWARD y:\n");CHKERRQ(ierr); 9897929ea7SJunchao Zhang ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 9997929ea7SJunchao Zhang 10097929ea7SJunchao Zhang /* Test reverse vecscatter */ 10197929ea7SJunchao Zhang ierr = VecSet(x,0.0);CHKERRQ(ierr); 10297929ea7SJunchao Zhang ierr = VecSet(y,0.0);CHKERRQ(ierr); 10397929ea7SJunchao Zhang ierr = VecGetOwnershipRange(y,&low,NULL);CHKERRQ(ierr); 10497929ea7SJunchao Zhang ierr = VecGetArray(y,&array);CHKERRQ(ierr); 10597929ea7SJunchao Zhang for (i=0; i<n; i++) { 10697929ea7SJunchao Zhang array[i] = (PetscScalar)(i+ low); 10797929ea7SJunchao Zhang } 10897929ea7SJunchao Zhang ierr = VecRestoreArray(y,&array);CHKERRQ(ierr); 10997929ea7SJunchao Zhang ierr = VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 11097929ea7SJunchao Zhang ierr = VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 11197929ea7SJunchao Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\nSCATTER_REVERSE x:\n");CHKERRQ(ierr); 11297929ea7SJunchao Zhang ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 11397929ea7SJunchao Zhang 11497929ea7SJunchao Zhang /* Free objects */ 11597929ea7SJunchao Zhang ierr = VecScatterDestroy(&ctx);CHKERRQ(ierr); 11697929ea7SJunchao Zhang ierr = ISDestroy(&isx);CHKERRQ(ierr); 11797929ea7SJunchao Zhang ierr = ISDestroy(&isy);CHKERRQ(ierr); 11897929ea7SJunchao Zhang ierr = VecDestroy(&x);CHKERRQ(ierr); 11997929ea7SJunchao Zhang ierr = VecDestroy(&x1);CHKERRQ(ierr); 12097929ea7SJunchao Zhang ierr = VecDestroy(&y);CHKERRQ(ierr); 12197929ea7SJunchao Zhang ierr = PetscFinalize(); 12297929ea7SJunchao Zhang return ierr; 12397929ea7SJunchao Zhang } 12497929ea7SJunchao Zhang 12597929ea7SJunchao Zhang /*TEST 12697929ea7SJunchao Zhang 12797929ea7SJunchao Zhang test: 12897929ea7SJunchao Zhang nsize: 2 12997929ea7SJunchao Zhang 13097929ea7SJunchao Zhang TEST*/ 131