xref: /petsc/src/vec/is/sf/tests/ex12.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
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   PetscInt       bs=1,n=5,i,low;
1097929ea7SJunchao Zhang   PetscInt       ix0[3] = {5,7,9},iy0[3] = {1,2,4},ix1[3] = {2,3,4},iy1[3] = {0,1,3};
1197929ea7SJunchao Zhang   PetscMPIInt    size,rank;
1297929ea7SJunchao Zhang   PetscScalar    *array;
1397929ea7SJunchao Zhang   Vec            x,y;
1497929ea7SJunchao Zhang   IS             isx,isy;
1597929ea7SJunchao Zhang   VecScatter     ctx;
1697929ea7SJunchao Zhang 
179566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
2097929ea7SJunchao Zhang 
21*08401ef6SPierre Jolivet   PetscCheck(size >=2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run more than one processor");
2297929ea7SJunchao Zhang 
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL));
2497929ea7SJunchao Zhang   n    = bs*n;
2597929ea7SJunchao Zhang 
2697929ea7SJunchao Zhang   /* Create vector x over shared memory */
279566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
289566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x,n,PETSC_DECIDE));
299566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
3097929ea7SJunchao Zhang 
319566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x,&low,NULL));
329566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x,&array));
3397929ea7SJunchao Zhang   for (i=0; i<n; i++) {
3497929ea7SJunchao Zhang     array[i] = (PetscScalar)(i + low);
3597929ea7SJunchao Zhang   }
369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x,&array));
3797929ea7SJunchao Zhang 
3897929ea7SJunchao Zhang   /* Create a sequential vector y */
399566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&y));
409566063dSJacob Faibussowitsch   PetscCall(VecSet(y,0.0));
4197929ea7SJunchao Zhang 
4297929ea7SJunchao Zhang   /* Create two index sets */
43dd400576SPatrick Sanan   if (rank == 0) {
449566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx));
459566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy));
4697929ea7SJunchao Zhang   } else {
479566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx));
489566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy));
4997929ea7SJunchao Zhang   }
5097929ea7SJunchao Zhang 
5197929ea7SJunchao Zhang   if (rank == 10) {
529566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank));
539566063dSJacob Faibussowitsch     PetscCall(ISView(isx,PETSC_VIEWER_STDOUT_SELF));
5497929ea7SJunchao Zhang   }
5597929ea7SJunchao Zhang 
569566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x,isx,y,isy,&ctx));
579566063dSJacob Faibussowitsch   PetscCall(VecScatterSetFromOptions(ctx));
5897929ea7SJunchao Zhang 
5997929ea7SJunchao Zhang   /* Test forward vecscatter */
609566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD));
619566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD));
6297929ea7SJunchao Zhang   if (rank == 0) {
639566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] y:\n",rank));
649566063dSJacob Faibussowitsch     PetscCall(VecView(y,PETSC_VIEWER_STDOUT_SELF));
6597929ea7SJunchao Zhang   }
6697929ea7SJunchao Zhang 
6797929ea7SJunchao Zhang   /* Test reverse vecscatter */
689566063dSJacob Faibussowitsch   PetscCall(VecScale(y,-1.0));
6997929ea7SJunchao Zhang   if (rank) {
709566063dSJacob Faibussowitsch     PetscCall(VecScale(y,1.0/(size - 1)));
7197929ea7SJunchao Zhang   }
7297929ea7SJunchao Zhang 
739566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_REVERSE));
749566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_REVERSE));
759566063dSJacob Faibussowitsch   PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
7697929ea7SJunchao Zhang 
7797929ea7SJunchao Zhang   /* Free spaces */
789566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&ctx));
799566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isx));
809566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isy));
819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
839566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
84b122ec5aSJacob Faibussowitsch   return 0;
8597929ea7SJunchao Zhang }
8697929ea7SJunchao Zhang 
8797929ea7SJunchao Zhang /*TEST
8897929ea7SJunchao Zhang 
8997929ea7SJunchao Zhang    test:
9097929ea7SJunchao Zhang       nsize: 3
9197929ea7SJunchao Zhang 
9297929ea7SJunchao Zhang TEST*/
93