xref: /petsc/src/vec/is/sf/tests/ex11.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
1 
2 static char help[]= "Scatters between parallel vectors. \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,N,i,low;
10   PetscInt       ix0[3] = {5,7,9},iy0[3] = {1,2,4},ix1[3] = {2,3,1},iy1[3] = {0,3,9};
11   PetscMPIInt    size,rank;
12   PetscScalar    *array;
13   Vec            x,x1,y;
14   IS             isx,isy;
15   VecScatter     ctx;
16   VecScatterType type;
17   PetscBool      flg;
18 
19   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
20   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
21   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
22 
23   PetscCheckFalse(size <2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run more than one processor");
24 
25   PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL));
26   n    = bs*n;
27 
28   /* Create vector x over shared memory */
29   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
30   PetscCall(VecSetSizes(x,n,PETSC_DECIDE));
31   PetscCall(VecSetFromOptions(x));
32 
33   PetscCall(VecGetOwnershipRange(x,&low,NULL));
34   PetscCall(VecGetArray(x,&array));
35   for (i=0; i<n; i++) {
36     array[i] = (PetscScalar)(i + low);
37   }
38   PetscCall(VecRestoreArray(x,&array));
39   /* PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); */
40 
41   /* Test some vector functions */
42   PetscCall(VecAssemblyBegin(x));
43   PetscCall(VecAssemblyEnd(x));
44 
45   PetscCall(VecGetSize(x,&N));
46   PetscCall(VecGetLocalSize(x,&n));
47 
48   PetscCall(VecDuplicate(x,&x1));
49   PetscCall(VecCopy(x,x1));
50   PetscCall(VecEqual(x,x1,&flg));
51   PetscCheck(flg,PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_WRONG,"x1 != x");
52 
53   PetscCall(VecScale(x1,2.0));
54   PetscCall(VecSet(x1,10.0));
55   /* PetscCall(VecView(x1,PETSC_VIEWER_STDOUT_WORLD)); */
56 
57   /* Create vector y over shared memory */
58   PetscCall(VecCreate(PETSC_COMM_WORLD,&y));
59   PetscCall(VecSetSizes(y,n,PETSC_DECIDE));
60   PetscCall(VecSetFromOptions(y));
61   PetscCall(VecGetArray(y,&array));
62   for (i=0; i<n; i++) {
63     array[i] = -(PetscScalar) (i + 100*rank);
64   }
65   PetscCall(VecRestoreArray(y,&array));
66   PetscCall(VecAssemblyBegin(y));
67   PetscCall(VecAssemblyEnd(y));
68   /* PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); */
69 
70   /* Create two index sets */
71   if (rank == 0) {
72     PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx));
73     PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy));
74   } else {
75     PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx));
76     PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy));
77   }
78 
79   if (rank == 10) {
80     PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank));
81     PetscCall(ISView(isx,PETSC_VIEWER_STDOUT_SELF));
82     PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n[%d] isy:\n",rank));
83     PetscCall(ISView(isy,PETSC_VIEWER_STDOUT_SELF));
84   }
85 
86   /* Create Vector scatter */
87   PetscCall(VecScatterCreate(x,isx,y,isy,&ctx));
88   PetscCall(VecScatterSetFromOptions(ctx));
89   PetscCall(VecScatterGetType(ctx,&type));
90   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"scatter type %s\n",type));
91 
92   /* Test forward vecscatter */
93   PetscCall(VecSet(y,0.0));
94   PetscCall(VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD));
95   PetscCall(VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD));
96   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSCATTER_FORWARD y:\n"));
97   PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
98 
99   /* Test reverse vecscatter */
100   PetscCall(VecSet(x,0.0));
101   PetscCall(VecSet(y,0.0));
102   PetscCall(VecGetOwnershipRange(y,&low,NULL));
103   PetscCall(VecGetArray(y,&array));
104   for (i=0; i<n; i++) {
105     array[i] = (PetscScalar)(i+ low);
106   }
107   PetscCall(VecRestoreArray(y,&array));
108   PetscCall(VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_REVERSE));
109   PetscCall(VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_REVERSE));
110   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSCATTER_REVERSE x:\n"));
111   PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
112 
113   /* Free objects */
114   PetscCall(VecScatterDestroy(&ctx));
115   PetscCall(ISDestroy(&isx));
116   PetscCall(ISDestroy(&isy));
117   PetscCall(VecDestroy(&x));
118   PetscCall(VecDestroy(&x1));
119   PetscCall(VecDestroy(&y));
120   PetscCall(PetscFinalize());
121   return 0;
122 }
123 
124 /*TEST
125 
126    test:
127       nsize: 2
128 
129 TEST*/
130