Actual source code: ex24.c
petsc-3.4.2 2013-07-02
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(NULL,"-n",&n,NULL);
23: PetscOptionsGetInt(NULL,"-bs",&bs,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) m = n + 2;
35: else m = n;
37: PetscMalloc((m)*sizeof(PetscInt),&blks);
38: blks[0] = n*rank;
39: for (i=1; i<m; i++) blks[i] = blks[i-1] + 1;
40: ISCreateBlock(PETSC_COMM_SELF,bs,m,blks,PETSC_COPY_VALUES,&is1);
41: PetscFree(blks);
43: VecCreateSeq(PETSC_COMM_SELF,bs*m,&y);
44: ISCreateStride(PETSC_COMM_SELF,bs*m,0,1,&is2);
46: /* each processor inserts the entire vector */
47: /* this is redundant but tests assembly */
48: for (i=0; i<bs*n*size; i++) {
49: value = (PetscScalar) i;
50: VecSetValues(x,1,&i,&value,INSERT_VALUES);
51: }
52: VecAssemblyBegin(x);
53: VecAssemblyEnd(x);
54: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
56: VecScatterCreate(x,is1,y,is2,&ctx);
57: VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);
58: VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);
60: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"----\n");
61: PetscViewerGetSingleton(PETSC_VIEWER_STDOUT_WORLD,&sviewer);
62: VecView(y,sviewer); fflush(stdout);
63: PetscViewerRestoreSingleton(PETSC_VIEWER_STDOUT_WORLD,&sviewer);
64: PetscSynchronizedFlush(PETSC_COMM_WORLD);
66: VecScatterDestroy(&ctx);
68: VecDestroy(&x);
69: VecDestroy(&y);
70: ISDestroy(&is1);
71: ISDestroy(&is2);
73: PetscFinalize();
74: return 0;
75: }