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