xref: /petsc/src/vec/is/sf/tests/ex12.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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;
19*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
20*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
2197929ea7SJunchao Zhang 
222c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size <2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run more than one processor");
2397929ea7SJunchao Zhang 
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL));
2597929ea7SJunchao Zhang   n    = bs*n;
2697929ea7SJunchao Zhang 
2797929ea7SJunchao Zhang   /* Create vector x over shared memory */
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x));
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(x,n,PETSC_DECIDE));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(x));
3197929ea7SJunchao Zhang 
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(x,&low,NULL));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(x,&array));
3497929ea7SJunchao Zhang   for (i=0; i<n; i++) {
3597929ea7SJunchao Zhang     array[i] = (PetscScalar)(i + low);
3697929ea7SJunchao Zhang   }
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(x,&array));
3897929ea7SJunchao Zhang 
3997929ea7SJunchao Zhang   /* Create a sequential vector y */
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&y));
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(y,0.0));
4297929ea7SJunchao Zhang 
4397929ea7SJunchao Zhang   /* Create two index sets */
44dd400576SPatrick Sanan   if (rank == 0) {
45*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx));
46*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy));
4797929ea7SJunchao Zhang   } else {
48*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx));
49*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy));
5097929ea7SJunchao Zhang   }
5197929ea7SJunchao Zhang 
5297929ea7SJunchao Zhang   if (rank == 10) {
53*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank));
54*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISView(isx,PETSC_VIEWER_STDOUT_SELF));
5597929ea7SJunchao Zhang   }
5697929ea7SJunchao Zhang 
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(x,isx,y,isy,&ctx));
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterSetFromOptions(ctx));
5997929ea7SJunchao Zhang 
6097929ea7SJunchao Zhang   /* Test forward vecscatter */
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD));
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD));
6397929ea7SJunchao Zhang   if (rank == 0) {
64*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[%d] y:\n",rank));
65*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_SELF));
6697929ea7SJunchao Zhang   }
6797929ea7SJunchao Zhang 
6897929ea7SJunchao Zhang   /* Test reverse vecscatter */
69*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(y,-1.0));
7097929ea7SJunchao Zhang   if (rank) {
71*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScale(y,1.0/(size - 1)));
7297929ea7SJunchao Zhang   }
7397929ea7SJunchao Zhang 
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_REVERSE));
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_REVERSE));
76*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
7797929ea7SJunchao Zhang 
7897929ea7SJunchao Zhang   /* Free spaces */
79*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&ctx));
80*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isx));
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isy));
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
83*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
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