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