Actual source code: ex10.c
2: static char help[]= "Scatters from a parallel vector to a sequential vector.\n\
3: uses block index sets\n\n";
5: #include petsc.h
6: #include petscis.h
7: #include petscvec.h
8: #include petscsys.h
12: int main(int argc,char **argv)
13: {
15: PetscInt bs = 1,n = 5,ix0[3] = {5,7,9},ix1[3] = {2,3,4},i,iy0[3] = {1,2,4},iy1[3] = {0,1,3};
16: PetscMPIInt size,rank;
17: PetscScalar value;
18: Vec x,y;
19: IS isx,isy;
20: VecScatter ctx = 0,newctx;
22: PetscInitialize(&argc,&argv,(char*)0,help);
23: MPI_Comm_size(PETSC_COMM_WORLD,&size);
24: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
26: if (size != 2) SETERRQ(1,"Must run with 2 processors");
28: PetscOptionsGetInt(PETSC_NULL,"-bs",&bs,PETSC_NULL);
29: n = bs*n;
31: /* create two vectors */
32: VecCreate(PETSC_COMM_WORLD,&x);
33: VecSetSizes(x,PETSC_DECIDE,size*n);
34: VecSetFromOptions(x);
35: VecCreateSeq(PETSC_COMM_SELF,n,&y);
37: /* create two index sets */
38: for (i=0; i<3; i++) {
39: ix0[i] *= bs; ix1[i] *= bs;
40: iy0[i] *= bs; iy1[i] *= bs;
41: }
43: if (!rank) {
44: ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,&isx);
45: ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,&isy);
46: } else {
47: ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,&isx);
48: ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,&isy);
49: }
51: /* fill local part of parallel vector */
52: for (i=n*rank; i<n*(rank+1); i++) {
53: value = (PetscScalar) i;
54: VecSetValues(x,1,&i,&value,INSERT_VALUES);
55: }
56: VecAssemblyBegin(x);
57: VecAssemblyEnd(x);
59: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
61: /* fill local part of parallel vector */
62: for (i=0; i<n; i++) {
63: value = -(PetscScalar) (i + 100*rank);
64: VecSetValues(y,1,&i,&value,INSERT_VALUES);
65: }
66: VecAssemblyBegin(y);
67: VecAssemblyEnd(y);
70: VecScatterCreate(x,isx,y,isy,&ctx);
71: VecScatterCopy(ctx,&newctx);
72: VecScatterDestroy(ctx);
74: VecScatterBegin(newctx,y,x,INSERT_VALUES,SCATTER_REVERSE);
75: VecScatterEnd(newctx,y,x,INSERT_VALUES,SCATTER_REVERSE);
76: VecScatterDestroy(newctx);
78: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
80: ISDestroy(isx);
81: ISDestroy(isy);
82: VecDestroy(x);
83: VecDestroy(y);
85: PetscFinalize();
86: return 0;
87: }
88: