xref: /petsc/src/vec/is/sf/tests/ex12.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 
2 static char help[]= "Scatters from a parallel vector to a sequential vector. \n\
3 uses block index sets\n\n";
4 
5 #include <petscvec.h>
6 
7 int main(int argc,char **argv)
8 {
9   PetscErrorCode ierr;
10   PetscInt       bs=1,n=5,i,low;
11   PetscInt       ix0[3] = {5,7,9},iy0[3] = {1,2,4},ix1[3] = {2,3,4},iy1[3] = {0,1,3};
12   PetscMPIInt    size,rank;
13   PetscScalar    *array;
14   Vec            x,y;
15   IS             isx,isy;
16   VecScatter     ctx;
17 
18   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
19   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
20   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
21 
22   PetscCheckFalse(size <2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run more than one processor");
23 
24   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL));
25   n    = bs*n;
26 
27   /* Create vector x over shared memory */
28   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x));
29   CHKERRQ(VecSetSizes(x,n,PETSC_DECIDE));
30   CHKERRQ(VecSetFromOptions(x));
31 
32   CHKERRQ(VecGetOwnershipRange(x,&low,NULL));
33   CHKERRQ(VecGetArray(x,&array));
34   for (i=0; i<n; i++) {
35     array[i] = (PetscScalar)(i + low);
36   }
37   CHKERRQ(VecRestoreArray(x,&array));
38 
39   /* Create a sequential vector y */
40   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&y));
41   CHKERRQ(VecSet(y,0.0));
42 
43   /* Create two index sets */
44   if (rank == 0) {
45     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx));
46     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy));
47   } else {
48     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx));
49     CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy));
50   }
51 
52   if (rank == 10) {
53     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank));
54     CHKERRQ(ISView(isx,PETSC_VIEWER_STDOUT_SELF));
55   }
56 
57   CHKERRQ(VecScatterCreate(x,isx,y,isy,&ctx));
58   CHKERRQ(VecScatterSetFromOptions(ctx));
59 
60   /* Test forward vecscatter */
61   CHKERRQ(VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD));
62   CHKERRQ(VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD));
63   if (rank == 0) {
64     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[%d] y:\n",rank));
65     CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_SELF));
66   }
67 
68   /* Test reverse vecscatter */
69   CHKERRQ(VecScale(y,-1.0));
70   if (rank) {
71     CHKERRQ(VecScale(y,1.0/(size - 1)));
72   }
73 
74   CHKERRQ(VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_REVERSE));
75   CHKERRQ(VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_REVERSE));
76   CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
77 
78   /* Free spaces */
79   CHKERRQ(VecScatterDestroy(&ctx));
80   CHKERRQ(ISDestroy(&isx));
81   CHKERRQ(ISDestroy(&isy));
82   CHKERRQ(VecDestroy(&x));
83   CHKERRQ(VecDestroy(&y));
84   ierr = PetscFinalize();
85   return ierr;
86 }
87 
88 /*TEST
89 
90    test:
91       nsize: 3
92 
93 TEST*/
94