Actual source code: ex10.c
1: /*$Id: ex10.c,v 1.17 2001/03/23 23:21:30 balay Exp $*/
3: static char help[]= "Scatters from a parallel vector to a sequential vector.n
4: uses block index setsnn";
6: #include petsc.h
7: #include petscis.h
8: #include petscvec.h
9: #include petscsys.h
11: int main(int argc,char **argv)
12: {
13: int bs = 1,n = 5,ierr,ix0[3] = {5,7,9},ix1[3] = {2,3,4};
14: int size,rank,i,iy0[3] = {1,2,4},iy1[3] = {0,1,3};
15: Scalar value;
16: Vec x,y;
17: IS isx,isy;
18: VecScatter ctx = 0,newctx;
20: PetscInitialize(&argc,&argv,(char*)0,help);
21: MPI_Comm_size(PETSC_COMM_WORLD,&size);
22: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
24: if (size != 2) SETERRQ(1,"Must run with 2 processors");
26: PetscOptionsGetInt(PETSC_NULL,"-bs",&bs,PETSC_NULL);
27: n = bs*n;
29: /* create two vectors */
30: VecCreateMPI(PETSC_COMM_WORLD,PETSC_DECIDE,size*n,&x);
31: VecCreateSeq(PETSC_COMM_SELF,n,&y);
33: /* create two index sets */
34: for (i=0; i<3; i++) {
35: ix0[i] *= bs; ix1[i] *= bs;
36: iy0[i] *= bs; iy1[i] *= bs;
37: }
39: if (!rank) {
40: ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,&isx);
41: ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,&isy);
42: } else {
43: ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,&isx);
44: ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,&isy);
45: }
47: /* fill local part of parallel vector */
48: for (i=n*rank; i<n*(rank+1); i++) {
49: value = (Scalar) i;
50: VecSetValues(x,1,&i,&value,INSERT_VALUES);
51: }
52: VecAssemblyBegin(x);
53: VecAssemblyEnd(x);
55: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
57: /* fill local part of parallel vector */
58: for (i=0; i<n; i++) {
59: value = -(Scalar) (i + 100*rank);
60: VecSetValues(y,1,&i,&value,INSERT_VALUES);
61: }
62: VecAssemblyBegin(y);
63: VecAssemblyEnd(y);
66: VecScatterCreate(x,isx,y,isy,&ctx);
67: VecScatterCopy(ctx,&newctx);
68: VecScatterDestroy(ctx);
70: VecScatterBegin(y,x,INSERT_VALUES,SCATTER_REVERSE,newctx);
71: VecScatterEnd(y,x,INSERT_VALUES,SCATTER_REVERSE,newctx);
72: VecScatterDestroy(newctx);
74: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
76: ISDestroy(isx);
77: ISDestroy(isy);
78: VecDestroy(x);
79: VecDestroy(y);
81: PetscFinalize();
82: return 0;
83: }
84: