Actual source code: ex8.c
1:
2: static char help[] = "Demonstrates generating a slice from a DA Vector.\n\n";
4: #include petscda.h
5: #include petscsys.h
6: #include petscao.h
10: /*
11: Given a DA generates a VecScatter context that will deliver a slice
12: of the global vector to each processor. In this example, each processor
13: receives the values i=*, j=*, k=rank, i.e. one z plane.
15: Note: This code is written assuming only one degree of freedom per node.
16: For multiple degrees of freedom per node use ISCreateBlock()
17: instead of ISCreateGeneral().
18: */
19: PetscErrorCode GenerateSliceScatter(DA da,VecScatter *scatter,Vec *vslice)
20: {
21: AO ao;
22: PetscInt M,N,P,nslice,*sliceindices,count,i,j;
23: PetscMPIInt rank;
25: MPI_Comm comm;
26: Vec vglobal;
27: IS isfrom,isto;
29: PetscObjectGetComm((PetscObject)da,&comm);
30: MPI_Comm_rank(comm,&rank);
32: DAGetAO(da,&ao);
33: DAGetInfo(da,0,&M,&N,&P,0,0,0,0,0,0,0);
35: /*
36: nslice is number of degrees of freedom in this processors slice
37: if there are more processors then z plans the extra processors get 0
38: elements in their slice.
39: */
40: if (rank < P) {nslice = M*N;} else nslice = 0;
42: /*
43: Generate the local vector to hold this processors slice
44: */
45: VecCreateSeq(PETSC_COMM_SELF,nslice,vslice);
46: DACreateGlobalVector(da,&vglobal);
48: /*
49: Generate the indices for the slice in the "natural" global ordering
50: Note: this is just an example, one could select any subset of nodes
51: on each processor. Just list them in the global natural ordering.
53: */
54: PetscMalloc((nslice+1)*sizeof(PetscInt),&sliceindices);
55: count = 0;
56: if (rank < P) {
57: for (j=0; j<N; j++) {
58: for (i=0; i<M; i++) {
59: sliceindices[count++] = rank*M*N + j*M + i;
60: }
61: }
62: }
63: /*
64: Convert the indices to the "PETSc" global ordering
65: */
66: AOApplicationToPetsc(ao,nslice,sliceindices);
67:
68: /* Create the "from" and "to" index set */
69: /* This is to scatter from the global vector */
70: ISCreateGeneral(PETSC_COMM_SELF,nslice,sliceindices,&isfrom);
71: /* This is to gather into the local vector */
72: ISCreateStride(PETSC_COMM_SELF,nslice,0,1,&isto);
73: VecScatterCreate(vglobal,isfrom,*vslice,isto,scatter);
74: ISDestroy(isfrom);
75: ISDestroy(isto);
76: PetscFree(sliceindices);
77: return 0;
78: }
83: int main(int argc,char **argv)
84: {
85: PetscMPIInt rank;
86: PetscInt m = PETSC_DECIDE,n = PETSC_DECIDE,p = PETSC_DECIDE,M = 3,N = 5,P=3,s=1;
87: PetscInt *lx = PETSC_NULL,*ly = PETSC_NULL,*lz = PETSC_NULL;
89: PetscTruth flg;
90: DA da;
91: Vec local,global,vslice;
92: PetscScalar value;
93: DAPeriodicType wrap = DA_XYPERIODIC;
94: DAStencilType stencil_type = DA_STENCIL_BOX;
95: VecScatter scatter;
97: PetscInitialize(&argc,&argv,(char*)0,help);
98: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
100: /* Read options */
101: PetscOptionsGetInt(PETSC_NULL,"-M",&M,PETSC_NULL);
102: PetscOptionsGetInt(PETSC_NULL,"-N",&N,PETSC_NULL);
103: PetscOptionsGetInt(PETSC_NULL,"-P",&P,PETSC_NULL);
104: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
105: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
106: PetscOptionsGetInt(PETSC_NULL,"-p",&p,PETSC_NULL);
107: PetscOptionsGetInt(PETSC_NULL,"-s",&s,PETSC_NULL);
108: PetscOptionsHasName(PETSC_NULL,"-star",&flg);
109: if (flg) stencil_type = DA_STENCIL_STAR;
111: /* Create distributed array and get vectors */
112: DACreate3d(PETSC_COMM_WORLD,wrap,stencil_type,M,N,P,m,n,p,1,s,
113: lx,ly,lz,&da);
114: DAView(da,PETSC_VIEWER_DRAW_WORLD);
115: DACreateGlobalVector(da,&global);
116: DACreateLocalVector(da,&local);
118: GenerateSliceScatter(da,&scatter,&vslice);
120: /* Put the value rank+1 into all locations of vslice and transfer back to global vector */
121: value = 1.0 + rank;
122: VecSet(vslice,value);
123: VecScatterBegin(scatter,vslice,global,INSERT_VALUES,SCATTER_REVERSE);
124: VecScatterEnd(scatter,vslice,global,INSERT_VALUES,SCATTER_REVERSE);
126: VecView(global,PETSC_VIEWER_DRAW_WORLD);
128: VecDestroy(local);
129: VecDestroy(global);
130: DADestroy(da);
131: PetscFinalize();
132: return 0;
133: }