Actual source code: ex24.c
2: static char help[] = "Scatters from a parallel vector to a sequential vector.\n\
3: Tests where the local part of the scatter is a copy.\n\n";
5: #include petscvec.h
9: int main(int argc,char **argv)
10: {
12: PetscMPIInt size,rank;
13: PetscInt n = 5,i,*blks,bs = 1,m = 2;
14: PetscScalar value;
15: Vec x,y;
16: IS is1,is2;
17: VecScatter ctx = 0;
18: PetscViewer sviewer;
20: PetscInitialize(&argc,&argv,(char*)0,help);
22: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
23: PetscOptionsGetInt(PETSC_NULL,"-bs",&bs,PETSC_NULL);
25: MPI_Comm_size(PETSC_COMM_WORLD,&size);
26: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
28: /* create two vectors */
29: VecCreate(PETSC_COMM_WORLD,&x);
30: VecSetSizes(x,PETSC_DECIDE,size*bs*n);
31: VecSetFromOptions(x);
33: /* create two index sets */
34: if (rank < size-1) {
35: m = n + 2;
36: } else {
37: m = n;
38: }
39: PetscMalloc((m)*sizeof(PetscInt),&blks);
40: blks[0] = n*rank*bs;
41: for (i=1; i<m; i++) {
42: blks[i] = blks[i-1] + bs;
43: }
44: ISCreateBlock(PETSC_COMM_SELF,bs,m,blks,&is1);
45: PetscFree(blks);
47: VecCreateSeq(PETSC_COMM_SELF,bs*m,&y);
48: ISCreateStride(PETSC_COMM_SELF,bs*m,0,1,&is2);
50: /* each processor inserts the entire vector */
51: /* this is redundant but tests assembly */
52: for (i=0; i<bs*n*size; i++) {
53: value = (PetscScalar) i;
54: VecSetValues(x,1,&i,&value,INSERT_VALUES);
55: }
56: VecAssemblyBegin(x);
57: VecAssemblyEnd(x);
58: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
60: VecScatterCreate(x,is1,y,is2,&ctx);
61: VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);
62: VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);
64: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"----\n");
65: PetscViewerGetSingleton(PETSC_VIEWER_STDOUT_WORLD,&sviewer);
66: VecView(y,sviewer); fflush(stdout);
67: PetscViewerRestoreSingleton(PETSC_VIEWER_STDOUT_WORLD,&sviewer);
68: PetscSynchronizedFlush(PETSC_COMM_WORLD);
70: VecScatterDestroy(ctx);
72: VecDestroy(x);
73: VecDestroy(y);
74: ISDestroy(is1);
75: ISDestroy(is2);
77: PetscFinalize();
78: return 0;
79: }
80: