1*97929ea7SJunchao Zhang 2*97929ea7SJunchao Zhang static char help[]= "Scatters between parallel vectors. \n\ 3*97929ea7SJunchao Zhang uses block index sets\n\n"; 4*97929ea7SJunchao Zhang 5*97929ea7SJunchao Zhang #include <petscvec.h> 6*97929ea7SJunchao Zhang 7*97929ea7SJunchao Zhang int main(int argc,char **argv) 8*97929ea7SJunchao Zhang { 9*97929ea7SJunchao Zhang PetscErrorCode ierr; 10*97929ea7SJunchao Zhang PetscInt bs=1,n=5,N,i,low; 11*97929ea7SJunchao Zhang PetscInt ix0[3] = {5,7,9},iy0[3] = {1,2,4},ix1[3] = {2,3,1},iy1[3] = {0,3,9}; 12*97929ea7SJunchao Zhang PetscMPIInt size,rank; 13*97929ea7SJunchao Zhang PetscScalar *array; 14*97929ea7SJunchao Zhang Vec x,x1,y; 15*97929ea7SJunchao Zhang IS isx,isy; 16*97929ea7SJunchao Zhang VecScatter ctx; 17*97929ea7SJunchao Zhang VecScatterType type; 18*97929ea7SJunchao Zhang PetscBool flg; 19*97929ea7SJunchao Zhang 20*97929ea7SJunchao Zhang ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 21*97929ea7SJunchao Zhang ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 22*97929ea7SJunchao Zhang ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 23*97929ea7SJunchao Zhang 24*97929ea7SJunchao Zhang if (size <2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run more than one processor"); 25*97929ea7SJunchao Zhang 26*97929ea7SJunchao Zhang ierr = PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);CHKERRQ(ierr); 27*97929ea7SJunchao Zhang n = bs*n; 28*97929ea7SJunchao Zhang 29*97929ea7SJunchao Zhang /* Create vector x over shared memory */ 30*97929ea7SJunchao Zhang ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 31*97929ea7SJunchao Zhang ierr = VecSetSizes(x,n,PETSC_DECIDE);CHKERRQ(ierr); 32*97929ea7SJunchao Zhang ierr = VecSetFromOptions(x);CHKERRQ(ierr); 33*97929ea7SJunchao Zhang 34*97929ea7SJunchao Zhang ierr = VecGetOwnershipRange(x,&low,NULL);CHKERRQ(ierr); 35*97929ea7SJunchao Zhang ierr = VecGetArray(x,&array);CHKERRQ(ierr); 36*97929ea7SJunchao Zhang for (i=0; i<n; i++) { 37*97929ea7SJunchao Zhang array[i] = (PetscScalar)(i + low); 38*97929ea7SJunchao Zhang } 39*97929ea7SJunchao Zhang ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 40*97929ea7SJunchao Zhang /* ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 41*97929ea7SJunchao Zhang 42*97929ea7SJunchao Zhang /* Test some vector functions */ 43*97929ea7SJunchao Zhang ierr = VecAssemblyBegin(x);CHKERRQ(ierr); 44*97929ea7SJunchao Zhang ierr = VecAssemblyEnd(x);CHKERRQ(ierr); 45*97929ea7SJunchao Zhang 46*97929ea7SJunchao Zhang ierr = VecGetSize(x,&N);CHKERRQ(ierr); 47*97929ea7SJunchao Zhang ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr); 48*97929ea7SJunchao Zhang 49*97929ea7SJunchao Zhang ierr = VecDuplicate(x,&x1);CHKERRQ(ierr); 50*97929ea7SJunchao Zhang ierr = VecCopy(x,x1);CHKERRQ(ierr); 51*97929ea7SJunchao Zhang ierr = VecEqual(x,x1,&flg);CHKERRQ(ierr); 52*97929ea7SJunchao Zhang if (!flg) SETERRQ(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_WRONG,"x1 != x"); 53*97929ea7SJunchao Zhang 54*97929ea7SJunchao Zhang ierr = VecScale(x1,2.0);CHKERRQ(ierr); 55*97929ea7SJunchao Zhang ierr = VecSet(x1,10.0);CHKERRQ(ierr); 56*97929ea7SJunchao Zhang /* ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 57*97929ea7SJunchao Zhang 58*97929ea7SJunchao Zhang /* Create vector y over shared memory */ 59*97929ea7SJunchao Zhang ierr = VecCreate(PETSC_COMM_WORLD,&y);CHKERRQ(ierr); 60*97929ea7SJunchao Zhang ierr = VecSetSizes(y,n,PETSC_DECIDE);CHKERRQ(ierr); 61*97929ea7SJunchao Zhang ierr = VecSetFromOptions(y);CHKERRQ(ierr); 62*97929ea7SJunchao Zhang ierr = VecGetArray(y,&array);CHKERRQ(ierr); 63*97929ea7SJunchao Zhang for (i=0; i<n; i++) { 64*97929ea7SJunchao Zhang array[i] = -(PetscScalar) (i + 100*rank); 65*97929ea7SJunchao Zhang } 66*97929ea7SJunchao Zhang ierr = VecRestoreArray(y,&array);CHKERRQ(ierr); 67*97929ea7SJunchao Zhang ierr = VecAssemblyBegin(y);CHKERRQ(ierr); 68*97929ea7SJunchao Zhang ierr = VecAssemblyEnd(y);CHKERRQ(ierr); 69*97929ea7SJunchao Zhang /* ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 70*97929ea7SJunchao Zhang 71*97929ea7SJunchao Zhang /* Create two index sets */ 72*97929ea7SJunchao Zhang if (!rank) { 73*97929ea7SJunchao Zhang ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr); 74*97929ea7SJunchao Zhang ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy);CHKERRQ(ierr); 75*97929ea7SJunchao Zhang } else { 76*97929ea7SJunchao Zhang ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr); 77*97929ea7SJunchao Zhang ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy);CHKERRQ(ierr); 78*97929ea7SJunchao Zhang } 79*97929ea7SJunchao Zhang 80*97929ea7SJunchao Zhang if (rank == 10) { 81*97929ea7SJunchao Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank);CHKERRQ(ierr); 82*97929ea7SJunchao Zhang ierr = ISView(isx,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 83*97929ea7SJunchao Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"\n[%d] isy:\n",rank);CHKERRQ(ierr); 84*97929ea7SJunchao Zhang ierr = ISView(isy,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 85*97929ea7SJunchao Zhang } 86*97929ea7SJunchao Zhang 87*97929ea7SJunchao Zhang /* Create Vector scatter */ 88*97929ea7SJunchao Zhang ierr = VecScatterCreate(x,isx,y,isy,&ctx);CHKERRQ(ierr); 89*97929ea7SJunchao Zhang ierr = VecScatterSetFromOptions(ctx);CHKERRQ(ierr); 90*97929ea7SJunchao Zhang ierr = VecScatterGetType(ctx,&type);CHKERRQ(ierr); 91*97929ea7SJunchao Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"scatter type %s\n",type);CHKERRQ(ierr); 92*97929ea7SJunchao Zhang 93*97929ea7SJunchao Zhang /* Test forward vecscatter */ 94*97929ea7SJunchao Zhang ierr = VecSet(y,0.0);CHKERRQ(ierr); 95*97929ea7SJunchao Zhang ierr = VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 96*97929ea7SJunchao Zhang ierr = VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 97*97929ea7SJunchao Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\nSCATTER_FORWARD y:\n");CHKERRQ(ierr); 98*97929ea7SJunchao Zhang ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 99*97929ea7SJunchao Zhang 100*97929ea7SJunchao Zhang /* Test reverse vecscatter */ 101*97929ea7SJunchao Zhang ierr = VecSet(x,0.0);CHKERRQ(ierr); 102*97929ea7SJunchao Zhang ierr = VecSet(y,0.0);CHKERRQ(ierr); 103*97929ea7SJunchao Zhang ierr = VecGetOwnershipRange(y,&low,NULL);CHKERRQ(ierr); 104*97929ea7SJunchao Zhang ierr = VecGetArray(y,&array);CHKERRQ(ierr); 105*97929ea7SJunchao Zhang for (i=0; i<n; i++) { 106*97929ea7SJunchao Zhang array[i] = (PetscScalar)(i+ low); 107*97929ea7SJunchao Zhang } 108*97929ea7SJunchao Zhang ierr = VecRestoreArray(y,&array);CHKERRQ(ierr); 109*97929ea7SJunchao Zhang ierr = VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 110*97929ea7SJunchao Zhang ierr = VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 111*97929ea7SJunchao Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\nSCATTER_REVERSE x:\n");CHKERRQ(ierr); 112*97929ea7SJunchao Zhang ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 113*97929ea7SJunchao Zhang 114*97929ea7SJunchao Zhang /* Free objects */ 115*97929ea7SJunchao Zhang ierr = VecScatterDestroy(&ctx);CHKERRQ(ierr); 116*97929ea7SJunchao Zhang ierr = ISDestroy(&isx);CHKERRQ(ierr); 117*97929ea7SJunchao Zhang ierr = ISDestroy(&isy);CHKERRQ(ierr); 118*97929ea7SJunchao Zhang ierr = VecDestroy(&x);CHKERRQ(ierr); 119*97929ea7SJunchao Zhang ierr = VecDestroy(&x1);CHKERRQ(ierr); 120*97929ea7SJunchao Zhang ierr = VecDestroy(&y);CHKERRQ(ierr); 121*97929ea7SJunchao Zhang ierr = PetscFinalize(); 122*97929ea7SJunchao Zhang return ierr; 123*97929ea7SJunchao Zhang } 124*97929ea7SJunchao Zhang 125*97929ea7SJunchao Zhang /*TEST 126*97929ea7SJunchao Zhang 127*97929ea7SJunchao Zhang test: 128*97929ea7SJunchao Zhang nsize: 2 129*97929ea7SJunchao Zhang 130*97929ea7SJunchao Zhang TEST*/ 131