Actual source code: ex9.c
2: static char help[]= "Scatters from a parallel vector to a sequential vector.\n\n";
4: #include petscvec.h
5: #include petscsys.h
9: int main(int argc,char **argv)
10: {
12: PetscInt n = 5,i,idx2[3] = {0,2,3},idx1[3] = {0,1,2};
13: PetscMPIInt size,rank;
14: PetscScalar mone = -1.0,value;
15: Vec x,y;
16: IS is1,is2;
17: VecScatter ctx = 0;
19: PetscInitialize(&argc,&argv,(char*)0,help);
20: MPI_Comm_size(PETSC_COMM_WORLD,&size);
21: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
23: /* create two vectors */
24: VecCreate(PETSC_COMM_WORLD,&x);
25: VecSetSizes(x,PETSC_DECIDE,size*n);
26: VecSetFromOptions(x);
27: VecCreateSeq(PETSC_COMM_SELF,n,&y);
29: /* create two index sets */
30: ISCreateGeneral(PETSC_COMM_SELF,3,idx1,&is1);
31: ISCreateGeneral(PETSC_COMM_SELF,3,idx2,&is2);
33: /* fill local part of parallel vector */
34: for (i=n*rank; i<n*(rank+1); i++) {
35: value = (PetscScalar) i;
36: VecSetValues(x,1,&i,&value,INSERT_VALUES);
37: }
38: VecAssemblyBegin(x);
39: VecAssemblyEnd(x);
41: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
43: VecSet(y,mone);
45: VecScatterCreate(x,is1,y,is2,&ctx);
46: VecScatterBegin(x,y,INSERT_VALUES,SCATTER_FORWARD,ctx);
47: VecScatterEnd(x,y,INSERT_VALUES,SCATTER_FORWARD,ctx);
48: VecScatterDestroy(ctx);
50: if (!rank) {
51: PetscPrintf(PETSC_COMM_SELF,"scattered vector\n");
52: VecView(y,PETSC_VIEWER_STDOUT_SELF);
53: }
54: ISDestroy(is1);
55: ISDestroy(is2);
56: VecDestroy(x);
57: VecDestroy(y);
59: PetscFinalize();
60: return 0;
61: }
62: