xref: /petsc/src/vec/is/sf/tests/ex11.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
205f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
215f80ce2aSJacob 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 
255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL));
2697929ea7SJunchao Zhang   n    = bs*n;
2797929ea7SJunchao Zhang 
2897929ea7SJunchao Zhang   /* Create vector x over shared memory */
295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(x,n,PETSC_DECIDE));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(x));
3297929ea7SJunchao Zhang 
335f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(x,&low,NULL));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(x,&array));
3597929ea7SJunchao Zhang   for (i=0; i<n; i++) {
3697929ea7SJunchao Zhang     array[i] = (PetscScalar)(i + low);
3797929ea7SJunchao Zhang   }
385f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(x,&array));
395f80ce2aSJacob Faibussowitsch   /* CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); */
4097929ea7SJunchao Zhang 
4197929ea7SJunchao Zhang   /* Test some vector functions */
425f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(x));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(x));
4497929ea7SJunchao Zhang 
455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(x,&N));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(x,&n));
4797929ea7SJunchao Zhang 
485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&x1));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(x,x1));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(VecEqual(x,x1,&flg));
5128b400f6SJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_WRONG,"x1 != x");
5297929ea7SJunchao Zhang 
535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(x1,2.0));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x1,10.0));
555f80ce2aSJacob Faibussowitsch   /* CHKERRQ(VecView(x1,PETSC_VIEWER_STDOUT_WORLD)); */
5697929ea7SJunchao Zhang 
5797929ea7SJunchao Zhang   /* Create vector y over shared memory */
585f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&y));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(y,n,PETSC_DECIDE));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(y));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(y,&array));
6297929ea7SJunchao Zhang   for (i=0; i<n; i++) {
6397929ea7SJunchao Zhang     array[i] = -(PetscScalar) (i + 100*rank);
6497929ea7SJunchao Zhang   }
655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(y,&array));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(y));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(y));
685f80ce2aSJacob Faibussowitsch   /* CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); */
6997929ea7SJunchao Zhang 
7097929ea7SJunchao Zhang   /* Create two index sets */
71dd400576SPatrick Sanan   if (rank == 0) {
725f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx));
735f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy));
7497929ea7SJunchao Zhang   } else {
755f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx));
765f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy));
7797929ea7SJunchao Zhang   }
7897929ea7SJunchao Zhang 
7997929ea7SJunchao Zhang   if (rank == 10) {
805f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank));
815f80ce2aSJacob Faibussowitsch     CHKERRQ(ISView(isx,PETSC_VIEWER_STDOUT_SELF));
825f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n[%d] isy:\n",rank));
835f80ce2aSJacob Faibussowitsch     CHKERRQ(ISView(isy,PETSC_VIEWER_STDOUT_SELF));
8497929ea7SJunchao Zhang   }
8597929ea7SJunchao Zhang 
8697929ea7SJunchao Zhang   /* Create Vector scatter */
875f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(x,isx,y,isy,&ctx));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterSetFromOptions(ctx));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterGetType(ctx,&type));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"scatter type %s\n",type));
9197929ea7SJunchao Zhang 
9297929ea7SJunchao Zhang   /* Test forward vecscatter */
935f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(y,0.0));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD));
955f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nSCATTER_FORWARD y:\n"));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
9897929ea7SJunchao Zhang 
9997929ea7SJunchao Zhang   /* Test reverse vecscatter */
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x,0.0));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(y,0.0));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(y,&low,NULL));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(y,&array));
10497929ea7SJunchao Zhang   for (i=0; i<n; i++) {
10597929ea7SJunchao Zhang     array[i] = (PetscScalar)(i+ low);
10697929ea7SJunchao Zhang   }
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(y,&array));
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_REVERSE));
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_REVERSE));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nSCATTER_REVERSE x:\n"));
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
11297929ea7SJunchao Zhang 
11397929ea7SJunchao Zhang   /* Free objects */
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&ctx));
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isx));
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isy));
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x1));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
120*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
121*b122ec5aSJacob Faibussowitsch   return 0;
12297929ea7SJunchao Zhang }
12397929ea7SJunchao Zhang 
12497929ea7SJunchao Zhang /*TEST
12597929ea7SJunchao Zhang 
12697929ea7SJunchao Zhang    test:
12797929ea7SJunchao Zhang       nsize: 2
12897929ea7SJunchao Zhang 
12997929ea7SJunchao Zhang TEST*/
130