xref: /petsc/src/vec/is/sf/tests/ex13.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
197929ea7SJunchao Zhang 
297929ea7SJunchao Zhang static char help[]= "Scatters from a sequential vector to a parallel 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   PetscViewer    sviewer;
1897929ea7SJunchao Zhang 
1997929ea7SJunchao Zhang   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
20*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
21*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
2297929ea7SJunchao Zhang 
232c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size <2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run more than one processor");
2497929ea7SJunchao Zhang 
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL));
2697929ea7SJunchao Zhang   n    = bs*n;
2797929ea7SJunchao Zhang 
2897929ea7SJunchao Zhang   /* Create vector x over shared memory */
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(x,n,PETSC_DECIDE));
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(x));
3297929ea7SJunchao Zhang 
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(x,&low,NULL));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(x,&array));
3597929ea7SJunchao Zhang   for (i=0; i<n; i++) {
3697929ea7SJunchao Zhang     array[i] = (PetscScalar)(i + low);
3797929ea7SJunchao Zhang   }
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(x,&array));
3997929ea7SJunchao Zhang 
4097929ea7SJunchao Zhang   /* Create a sequential vector y */
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&y));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(y,&array));
4397929ea7SJunchao Zhang   for (i=0; i<n; i++) {
4497929ea7SJunchao Zhang     array[i] = (PetscScalar)(i + 100*rank);
4597929ea7SJunchao Zhang   }
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(y,&array));
4797929ea7SJunchao Zhang 
4897929ea7SJunchao Zhang   /* Create two index sets */
49dd400576SPatrick Sanan   if (rank == 0) {
50*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx));
51*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy));
5297929ea7SJunchao Zhang   } else {
53*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx));
54*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy));
5597929ea7SJunchao Zhang   }
5697929ea7SJunchao Zhang 
5797929ea7SJunchao Zhang   if (rank == 10) {
58*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank));
59*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISView(isx,PETSC_VIEWER_STDOUT_SELF));
6097929ea7SJunchao Zhang   }
6197929ea7SJunchao Zhang 
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(y,isy,x,isx,&ctx));
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterSetFromOptions(ctx));
6497929ea7SJunchao Zhang 
6597929ea7SJunchao Zhang   /* Test forward vecscatter */
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x,0.0));
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_FORWARD));
68*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_FORWARD));
69*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
7097929ea7SJunchao Zhang 
7197929ea7SJunchao Zhang   /* Test reverse vecscatter */
72*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(x,-1.0));
73*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_REVERSE));
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_REVERSE));
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer));
7697929ea7SJunchao Zhang   if (rank == 1) {
77*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(y,sviewer));
7897929ea7SJunchao Zhang   }
79*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer));
8097929ea7SJunchao Zhang 
8197929ea7SJunchao Zhang   /* Free spaces */
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&ctx));
83*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isx));
84*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isy));
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
86*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
8797929ea7SJunchao Zhang   ierr = PetscFinalize();
8897929ea7SJunchao Zhang   return ierr;
8997929ea7SJunchao Zhang }
9097929ea7SJunchao Zhang 
9197929ea7SJunchao Zhang /*TEST
9297929ea7SJunchao Zhang 
9397929ea7SJunchao Zhang    test:
9497929ea7SJunchao Zhang       nsize: 3
9597929ea7SJunchao Zhang 
9697929ea7SJunchao Zhang TEST*/
97