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