Actual source code: ex10.c
1: /*$Id: ex10.c,v 1.20 2001/09/11 16:32:10 bsmith 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: PetscScalar 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: VecCreate(PETSC_COMM_WORLD,&x);
31: VecSetSizes(x,PETSC_DECIDE,size*n);
32: VecSetFromOptions(x);
33: VecCreateSeq(PETSC_COMM_SELF,n,&y);
35: /* create two index sets */
36: for (i=0; i<3; i++) {
37: ix0[i] *= bs; ix1[i] *= bs;
38: iy0[i] *= bs; iy1[i] *= bs;
39: }
41: if (!rank) {
42: ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,&isx);
43: ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,&isy);
44: } else {
45: ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,&isx);
46: ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,&isy);
47: }
49: /* fill local part of parallel vector */
50: for (i=n*rank; i<n*(rank+1); i++) {
51: value = (PetscScalar) i;
52: VecSetValues(x,1,&i,&value,INSERT_VALUES);
53: }
54: VecAssemblyBegin(x);
55: VecAssemblyEnd(x);
57: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
59: /* fill local part of parallel vector */
60: for (i=0; i<n; i++) {
61: value = -(PetscScalar) (i + 100*rank);
62: VecSetValues(y,1,&i,&value,INSERT_VALUES);
63: }
64: VecAssemblyBegin(y);
65: VecAssemblyEnd(y);
68: VecScatterCreate(x,isx,y,isy,&ctx);
69: VecScatterCopy(ctx,&newctx);
70: VecScatterDestroy(ctx);
72: VecScatterBegin(y,x,INSERT_VALUES,SCATTER_REVERSE,newctx);
73: VecScatterEnd(y,x,INSERT_VALUES,SCATTER_REVERSE,newctx);
74: VecScatterDestroy(newctx);
76: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
78: ISDestroy(isx);
79: ISDestroy(isy);
80: VecDestroy(x);
81: VecDestroy(y);
83: PetscFinalize();
84: return 0;
85: }
86: