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