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 petscvec.h
9: int main(int argc,char **argv)
10: {
12: 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};
13: PetscMPIInt size,rank;
14: PetscScalar value;
15: Vec x,y;
16: IS isx,isy;
17: VecScatter ctx = 0,newctx;
19: PetscInitialize(&argc,&argv,(char*)0,help);
20: MPI_Comm_size(PETSC_COMM_WORLD,&size);
21: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
23: if (size != 2) SETERRQ(1,"Must run with 2 processors");
25: PetscOptionsGetInt(PETSC_NULL,"-bs",&bs,PETSC_NULL);
26: n = bs*n;
28: /* create two vectors */
29: VecCreate(PETSC_COMM_WORLD,&x);
30: VecSetSizes(x,PETSC_DECIDE,size*n);
31: VecSetFromOptions(x);
32: VecCreateSeq(PETSC_COMM_SELF,n,&y);
34: /* create two index sets */
35: for (i=0; i<3; i++) {
36: ix0[i] *= bs; ix1[i] *= bs;
37: iy0[i] *= bs; iy1[i] *= bs;
38: }
40: if (!rank) {
41: ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,&isx);
42: ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,&isy);
43: } else {
44: ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,&isx);
45: ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,&isy);
46: }
48: /* fill local part of parallel vector */
49: for (i=n*rank; i<n*(rank+1); i++) {
50: value = (PetscScalar) i;
51: VecSetValues(x,1,&i,&value,INSERT_VALUES);
52: }
53: VecAssemblyBegin(x);
54: VecAssemblyEnd(x);
56: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
58: /* fill local part of parallel vector */
59: for (i=0; i<n; i++) {
60: value = -(PetscScalar) (i + 100*rank);
61: VecSetValues(y,1,&i,&value,INSERT_VALUES);
62: }
63: VecAssemblyBegin(y);
64: VecAssemblyEnd(y);
67: VecScatterCreate(x,isx,y,isy,&ctx);
68: VecScatterCopy(ctx,&newctx);
69: VecScatterDestroy(ctx);
71: VecScatterBegin(newctx,y,x,INSERT_VALUES,SCATTER_REVERSE);
72: VecScatterEnd(newctx,y,x,INSERT_VALUES,SCATTER_REVERSE);
73: VecScatterDestroy(newctx);
75: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
77: ISDestroy(isx);
78: ISDestroy(isy);
79: VecDestroy(x);
80: VecDestroy(y);
82: PetscFinalize();
83: return 0;
84: }
85: