Actual source code: ex24.c
1: /*$Id: ex24.c,v 1.19 2001/04/10 19:35:02 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: Scalar 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: VecCreateMPI(PETSC_COMM_WORLD,PETSC_DECIDE,size*bs*n,&x);
29: /* create two index sets */
30: if (rank < size-1) {
31: m = n + 2;
32: } else {
33: m = n;
34: }
35: PetscMalloc((m)*sizeof(int),&blks);
36: blks[0] = n*rank*bs;
37: for (i=1; i<m; i++) {
38: blks[i] = blks[i-1] + bs;
39: }
40: ISCreateBlock(PETSC_COMM_SELF,bs,m,blks,&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 = (Scalar) 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(x,y,INSERT_VALUES,SCATTER_FORWARD,ctx);
58: VecScatterEnd(x,y,INSERT_VALUES,SCATTER_FORWARD,ctx);
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: }
76: