Actual source code: ex26.c
1: /*$Id: ex26.c,v 1.10 2001/03/23 23:21:30 balay Exp $*/
2: /*
4: Test program follows. Writing it I realised that
5: 1/ I built the pipeline object around an MPI-to-MPI vector scatter.
6: 2/ That necessitates the 'little trick' below.
7: 3/ This trick will always be necessary, since the code inside the
8: pipe is a critical section.
9: 4/ Hence I really should have used an MPI-to-Seq scatter.
10: 5/ This shouldn't be too hard to fix in the implementation you
11: guys are making,right? :-) <-- smiley just in case.
13: If this is not clear, I'll try to elaborate.
15: */
16: /* Example of pipeline code:
17: accumulation of the sum $s_p=sum_{qleq p} (q+1)^2$.
18: E.g., processor 3 computes 1^2+2^2+3^+4^2 = 30.
19: Every processor computes its term, then passes it on to the next.
20: */
21: #include petscvec.h
23: int main(int Argc,char **Args)
24: {
25: Vec src_v,tar_v,loc_v;
26: IS src_idx,tar_idx;
27: VecPipeline pipe;
28: MPI_Comm comm;
29: int size,rank,src_loc,tar_loc,ierr,zero_loc=0;
30: Scalar zero=0,my_value,*vec_values,*loc_ar;
32: PetscInitialize(&Argc,&Args,PETSC_NULL,PETSC_NULL);
34: comm = MPI_COMM_WORLD;
35: MPI_Comm_size(comm,&size);
36: MPI_Comm_rank(comm,&rank);
37:
38: /* Create the necessary vectors; one element per processor */
39: VecCreateMPI(comm,1,PETSC_DECIDE,&tar_v);
40: VecSet(&zero,tar_v);
41: VecCreateMPI(comm,1,PETSC_DECIDE,&src_v);
42: VecCreateSeq(MPI_COMM_SELF,1,&loc_v);
43: /* -- little trick: we need a distributed and a local vector
44: that share each other's data; see below for application */
45: VecGetArray(loc_v,&loc_ar);
46: VecPlaceArray(src_v,loc_ar);
48: /* Create the pipeline data: we write into our own location,
49: and read one location from the left */
50: tar_loc = rank;
51: if (tar_loc>0) src_loc = tar_loc-1; else src_loc = tar_loc;
52: ISCreateGeneral(MPI_COMM_SELF,1,&tar_loc,&tar_idx);
53: ISCreateGeneral(MPI_COMM_SELF,1,&src_loc,&src_idx);
54: VecPipelineCreate(comm,src_v,src_idx,tar_v,tar_idx,&pipe);
55: VecPipelineSetType(pipe,PIPELINE_SEQUENTIAL,PETSC_NULL);
56: VecPipelineSetup(pipe);
58: /* The actual pipe:
59: receive accumulated value from previous processor,
60: add the square of your own value, and send on. */
61: VecPipelineBegin(src_v,tar_v,INSERT_VALUES,SCATTER_FORWARD,PIPELINE_UP,pipe);
62: VecGetArray(tar_v,&vec_values);
63: my_value = vec_values[0] + (double)((rank+1)*(rank+1));
65: VecRestoreArray(tar_v,&vec_values);CHKERRQ(ierr)
66: /* -- little trick: we have to be able to call VecAssembly,
67: but since this code executed sequentially (critical section!),
68: we have a local vector with data aliased to the distributed one */
69: VecSetValues(loc_v,1,&zero_loc,&my_value,INSERT_VALUES);
70: VecAssemblyBegin(loc_v);
71: VecAssemblyEnd(loc_v);
72: VecPipelineEnd(src_v,tar_v,INSERT_VALUES,SCATTER_FORWARD,PIPELINE_UP,pipe);
74: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] value=%dn",rank,(int)PetscRealPart(my_value));
75: PetscSynchronizedFlush(PETSC_COMM_WORLD);
77: /* Clean up */
78: VecPipelineDestroy(pipe);
79: VecDestroy(src_v);
80: VecDestroy(tar_v);
81: VecDestroy(loc_v);
82: ISDestroy(src_idx);
83: ISDestroy(tar_idx);
85: PetscFinalize();
87: return 0;
88: }