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