xref: /petsc/src/vec/is/sf/tests/ex12.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
197929ea7SJunchao Zhang 
297929ea7SJunchao Zhang static char help[]= "Scatters from a parallel vector to a sequential vector. \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,i,low;
1197929ea7SJunchao Zhang   PetscInt       ix0[3] = {5,7,9},iy0[3] = {1,2,4},ix1[3] = {2,3,4},iy1[3] = {0,1,3};
1297929ea7SJunchao Zhang   PetscMPIInt    size,rank;
1397929ea7SJunchao Zhang   PetscScalar    *array;
1497929ea7SJunchao Zhang   Vec            x,y;
1597929ea7SJunchao Zhang   IS             isx,isy;
1697929ea7SJunchao Zhang   VecScatter     ctx;
1797929ea7SJunchao Zhang 
1897929ea7SJunchao Zhang   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
19ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
20ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
2197929ea7SJunchao Zhang 
22*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size <2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run more than one processor");
2397929ea7SJunchao Zhang 
2497929ea7SJunchao Zhang   ierr = PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);CHKERRQ(ierr);
2597929ea7SJunchao Zhang   n    = bs*n;
2697929ea7SJunchao Zhang 
2797929ea7SJunchao Zhang   /* Create vector x over shared memory */
2897929ea7SJunchao Zhang   ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
2997929ea7SJunchao Zhang   ierr = VecSetSizes(x,n,PETSC_DECIDE);CHKERRQ(ierr);
3097929ea7SJunchao Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
3197929ea7SJunchao Zhang 
3297929ea7SJunchao Zhang   ierr = VecGetOwnershipRange(x,&low,NULL);CHKERRQ(ierr);
3397929ea7SJunchao Zhang   ierr = VecGetArray(x,&array);CHKERRQ(ierr);
3497929ea7SJunchao Zhang   for (i=0; i<n; i++) {
3597929ea7SJunchao Zhang     array[i] = (PetscScalar)(i + low);
3697929ea7SJunchao Zhang   }
3797929ea7SJunchao Zhang   ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
3897929ea7SJunchao Zhang 
3997929ea7SJunchao Zhang   /* Create a sequential vector y */
4097929ea7SJunchao Zhang   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&y);CHKERRQ(ierr);
4197929ea7SJunchao Zhang   ierr = VecSet(y,0.0);CHKERRQ(ierr);
4297929ea7SJunchao Zhang 
4397929ea7SJunchao Zhang   /* Create two index sets */
44dd400576SPatrick Sanan   if (rank == 0) {
4597929ea7SJunchao Zhang     ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr);
4697929ea7SJunchao Zhang     ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy);CHKERRQ(ierr);
4797929ea7SJunchao Zhang   } else {
4897929ea7SJunchao Zhang     ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr);
4997929ea7SJunchao Zhang     ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy);CHKERRQ(ierr);
5097929ea7SJunchao Zhang   }
5197929ea7SJunchao Zhang 
5297929ea7SJunchao Zhang   if (rank == 10) {
5397929ea7SJunchao Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank);CHKERRQ(ierr);
5497929ea7SJunchao Zhang     ierr = ISView(isx,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
5597929ea7SJunchao Zhang   }
5697929ea7SJunchao Zhang 
5797929ea7SJunchao Zhang   ierr = VecScatterCreate(x,isx,y,isy,&ctx);CHKERRQ(ierr);
5897929ea7SJunchao Zhang   ierr = VecScatterSetFromOptions(ctx);CHKERRQ(ierr);
5997929ea7SJunchao Zhang 
6097929ea7SJunchao Zhang   /* Test forward vecscatter */
6197929ea7SJunchao Zhang   ierr = VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6297929ea7SJunchao Zhang   ierr = VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6397929ea7SJunchao Zhang   if (rank == 0) {
6497929ea7SJunchao Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] y:\n",rank);CHKERRQ(ierr);
6597929ea7SJunchao Zhang     ierr = VecView(y,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
6697929ea7SJunchao Zhang   }
6797929ea7SJunchao Zhang 
6897929ea7SJunchao Zhang   /* Test reverse vecscatter */
6997929ea7SJunchao Zhang   ierr = VecScale(y,-1.0);CHKERRQ(ierr);
7097929ea7SJunchao Zhang   if (rank) {
7197929ea7SJunchao Zhang     ierr = VecScale(y,1.0/(size - 1));CHKERRQ(ierr);
7297929ea7SJunchao Zhang   }
7397929ea7SJunchao Zhang 
7497929ea7SJunchao Zhang   ierr = VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7597929ea7SJunchao Zhang   ierr = VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7697929ea7SJunchao Zhang   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
7797929ea7SJunchao Zhang 
7897929ea7SJunchao Zhang   /* Free spaces */
7997929ea7SJunchao Zhang   ierr = VecScatterDestroy(&ctx);CHKERRQ(ierr);
8097929ea7SJunchao Zhang   ierr = ISDestroy(&isx);CHKERRQ(ierr);
8197929ea7SJunchao Zhang   ierr = ISDestroy(&isy);CHKERRQ(ierr);
8297929ea7SJunchao Zhang   ierr = VecDestroy(&x);CHKERRQ(ierr);
8397929ea7SJunchao Zhang   ierr = VecDestroy(&y);CHKERRQ(ierr);
8497929ea7SJunchao Zhang   ierr = PetscFinalize();
8597929ea7SJunchao Zhang   return ierr;
8697929ea7SJunchao Zhang }
8797929ea7SJunchao Zhang 
8897929ea7SJunchao Zhang /*TEST
8997929ea7SJunchao Zhang 
9097929ea7SJunchao Zhang    test:
9197929ea7SJunchao Zhang       nsize: 3
9297929ea7SJunchao Zhang 
9397929ea7SJunchao Zhang TEST*/
94