xref: /petsc/src/vec/is/sf/tests/ex11.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
197929ea7SJunchao Zhang 
297929ea7SJunchao Zhang static char help[]= "Scatters between parallel vectors. \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,N,i,low;
1097929ea7SJunchao Zhang   PetscInt       ix0[3] = {5,7,9},iy0[3] = {1,2,4},ix1[3] = {2,3,1},iy1[3] = {0,3,9};
1197929ea7SJunchao Zhang   PetscMPIInt    size,rank;
1297929ea7SJunchao Zhang   PetscScalar    *array;
1397929ea7SJunchao Zhang   Vec            x,x1,y;
1497929ea7SJunchao Zhang   IS             isx,isy;
1597929ea7SJunchao Zhang   VecScatter     ctx;
1697929ea7SJunchao Zhang   VecScatterType type;
1797929ea7SJunchao Zhang   PetscBool      flg;
1897929ea7SJunchao Zhang 
19*327415f7SBarry Smith   PetscFunctionBeginUser;
209566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
2397929ea7SJunchao Zhang 
2408401ef6SPierre Jolivet   PetscCheck(size >=2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run more than one processor");
2597929ea7SJunchao Zhang 
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL));
2797929ea7SJunchao Zhang   n    = bs*n;
2897929ea7SJunchao Zhang 
2997929ea7SJunchao Zhang   /* Create vector x over shared memory */
309566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
319566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x,n,PETSC_DECIDE));
329566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
3397929ea7SJunchao Zhang 
349566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x,&low,NULL));
359566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x,&array));
3697929ea7SJunchao Zhang   for (i=0; i<n; i++) {
3797929ea7SJunchao Zhang     array[i] = (PetscScalar)(i + low);
3897929ea7SJunchao Zhang   }
399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x,&array));
409566063dSJacob Faibussowitsch   /* PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); */
4197929ea7SJunchao Zhang 
4297929ea7SJunchao Zhang   /* Test some vector functions */
439566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(x));
449566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(x));
4597929ea7SJunchao Zhang 
469566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x,&N));
479566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(x,&n));
4897929ea7SJunchao Zhang 
499566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&x1));
509566063dSJacob Faibussowitsch   PetscCall(VecCopy(x,x1));
519566063dSJacob Faibussowitsch   PetscCall(VecEqual(x,x1,&flg));
5228b400f6SJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_WRONG,"x1 != x");
5397929ea7SJunchao Zhang 
549566063dSJacob Faibussowitsch   PetscCall(VecScale(x1,2.0));
559566063dSJacob Faibussowitsch   PetscCall(VecSet(x1,10.0));
569566063dSJacob Faibussowitsch   /* PetscCall(VecView(x1,PETSC_VIEWER_STDOUT_WORLD)); */
5797929ea7SJunchao Zhang 
5897929ea7SJunchao Zhang   /* Create vector y over shared memory */
599566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&y));
609566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(y,n,PETSC_DECIDE));
619566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(y));
629566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y,&array));
6397929ea7SJunchao Zhang   for (i=0; i<n; i++) {
6497929ea7SJunchao Zhang     array[i] = -(PetscScalar) (i + 100*rank);
6597929ea7SJunchao Zhang   }
669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y,&array));
679566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(y));
689566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(y));
699566063dSJacob Faibussowitsch   /* PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); */
7097929ea7SJunchao Zhang 
7197929ea7SJunchao Zhang   /* Create two index sets */
72dd400576SPatrick Sanan   if (rank == 0) {
739566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx));
749566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy));
7597929ea7SJunchao Zhang   } else {
769566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx));
779566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy));
7897929ea7SJunchao Zhang   }
7997929ea7SJunchao Zhang 
8097929ea7SJunchao Zhang   if (rank == 10) {
819566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank));
829566063dSJacob Faibussowitsch     PetscCall(ISView(isx,PETSC_VIEWER_STDOUT_SELF));
839566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n[%d] isy:\n",rank));
849566063dSJacob Faibussowitsch     PetscCall(ISView(isy,PETSC_VIEWER_STDOUT_SELF));
8597929ea7SJunchao Zhang   }
8697929ea7SJunchao Zhang 
8797929ea7SJunchao Zhang   /* Create Vector scatter */
889566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x,isx,y,isy,&ctx));
899566063dSJacob Faibussowitsch   PetscCall(VecScatterSetFromOptions(ctx));
909566063dSJacob Faibussowitsch   PetscCall(VecScatterGetType(ctx,&type));
919566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"scatter type %s\n",type));
9297929ea7SJunchao Zhang 
9397929ea7SJunchao Zhang   /* Test forward vecscatter */
949566063dSJacob Faibussowitsch   PetscCall(VecSet(y,0.0));
959566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD));
969566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD));
979566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSCATTER_FORWARD y:\n"));
989566063dSJacob Faibussowitsch   PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
9997929ea7SJunchao Zhang 
10097929ea7SJunchao Zhang   /* Test reverse vecscatter */
1019566063dSJacob Faibussowitsch   PetscCall(VecSet(x,0.0));
1029566063dSJacob Faibussowitsch   PetscCall(VecSet(y,0.0));
1039566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(y,&low,NULL));
1049566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y,&array));
10597929ea7SJunchao Zhang   for (i=0; i<n; i++) {
10697929ea7SJunchao Zhang     array[i] = (PetscScalar)(i+ low);
10797929ea7SJunchao Zhang   }
1089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y,&array));
1099566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_REVERSE));
1109566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_REVERSE));
1119566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSCATTER_REVERSE x:\n"));
1129566063dSJacob Faibussowitsch   PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
11397929ea7SJunchao Zhang 
11497929ea7SJunchao Zhang   /* Free objects */
1159566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&ctx));
1169566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isx));
1179566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isy));
1189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x1));
1209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
1219566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
122b122ec5aSJacob Faibussowitsch   return 0;
12397929ea7SJunchao Zhang }
12497929ea7SJunchao Zhang 
12597929ea7SJunchao Zhang /*TEST
12697929ea7SJunchao Zhang 
12797929ea7SJunchao Zhang    test:
12897929ea7SJunchao Zhang       nsize: 2
12997929ea7SJunchao Zhang 
13097929ea7SJunchao Zhang TEST*/
131