Actual source code: ex8.c

  1: /*$Id: ex8.c,v 1.28 2001/08/07 03:04:42 balay Exp $*/
  2: 
  3: static char help[] = "Demonstrates generating a slice from a DA Vector.nn";

 5:  #include petscda.h
 6:  #include petscsys.h
 7:  #include petscao.h

  9: /*
 10:     Given a DA generates a VecScatter context that will deliver a slice
 11:   of the global vector to each processor. In this example, each processor
 12:   receives the values i=*, j=*, k=rank, i.e. one z plane.

 14:   Note: This code is written assuming only one degree of freedom per node.
 15:   For multiple degrees of freedom per node use ISCreateBlock()
 16:   instead of ISCreateGeneral().
 17: */
 18: int GenerateSliceScatter(DA da,VecScatter *scatter,Vec *vslice)
 19: {
 20:   AO       ao;
 21:   int      M,N,P,nslice,rank,*sliceindices,count,ierr,i,j;
 22:   MPI_Comm comm;
 23:   Vec      vglobal;
 24:   IS       isfrom,isto;

 26:   PetscObjectGetComm((PetscObject)da,&comm);
 27:   MPI_Comm_rank(comm,&rank);

 29:   DAGetAO(da,&ao);
 30:   DAGetInfo(da,0,&M,&N,&P,0,0,0,0,0,0,0);

 32:   /* 
 33:      nslice is number of degrees of freedom in this processors slice
 34:    if there are more processors then z plans the extra processors get 0
 35:    elements in their slice.
 36:   */
 37:   if (rank < P) {nslice = M*N;} else nslice = 0;

 39:   /* 
 40:      Generate the local vector to hold this processors slice
 41:   */
 42:   VecCreateSeq(PETSC_COMM_SELF,nslice,vslice);
 43:   DACreateGlobalVector(da,&vglobal);

 45:   /*
 46:        Generate the indices for the slice in the "natural" global ordering
 47:      Note: this is just an example, one could select any subset of nodes 
 48:     on each processor. Just list them in the global natural ordering.

 50:   */
 51:   PetscMalloc((nslice+1)*sizeof(int),&sliceindices);
 52:   count = 0;
 53:   if (rank < P) {
 54:     for (j=0; j<N; j++) {
 55:       for (i=0; i<M; i++) {
 56:          sliceindices[count++] = rank*M*N + j*M + i;
 57:       }
 58:     }
 59:   }
 60:   /*
 61:       Convert the indices to the "PETSc" global ordering
 62:   */
 63:   AOApplicationToPetsc(ao,nslice,sliceindices);
 64: 
 65:   /* Create the "from" and "to" index set */
 66:   /* This is to scatter from the global vector */
 67:   ISCreateGeneral(PETSC_COMM_SELF,nslice,sliceindices,&isfrom);
 68:   /* This is to gather into the local vector */
 69:   ISCreateStride(PETSC_COMM_SELF,nslice,0,1,&isto);

 71:   VecScatterCreate(vglobal,isfrom,*vslice,isto,scatter);

 73:   ISDestroy(isfrom);
 74:   ISDestroy(isto);

 76:   PetscFree(sliceindices);
 77:   return 0;
 78: }


 81: int main(int argc,char **argv)
 82: {
 83:   int            rank,M = 3,N = 5,P=3,s=1;
 84:   int            m = PETSC_DECIDE,n = PETSC_DECIDE,p = PETSC_DECIDE,ierr;
 85:   int            *lx = PETSC_NULL,*ly = PETSC_NULL,*lz = PETSC_NULL;
 86:   PetscTruth     flg;
 87:   DA             da;
 88:   Vec            local,global,vslice;
 89:   PetscScalar    value;
 90:   DAPeriodicType wrap = DA_XYPERIODIC;
 91:   DAStencilType  stencil_type = DA_STENCIL_BOX;
 92:   VecScatter     scatter;

 94:   PetscInitialize(&argc,&argv,(char*)0,help);
 95:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

 97:   /* Read options */
 98:   PetscOptionsGetInt(PETSC_NULL,"-M",&M,PETSC_NULL);
 99:   PetscOptionsGetInt(PETSC_NULL,"-N",&N,PETSC_NULL);
100:   PetscOptionsGetInt(PETSC_NULL,"-P",&P,PETSC_NULL);
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,"-s",&s,PETSC_NULL);
105:   PetscOptionsHasName(PETSC_NULL,"-star",&flg);
106:   if (flg) stencil_type =  DA_STENCIL_STAR;

108:   /* Create distributed array and get vectors */
109:   DACreate3d(PETSC_COMM_WORLD,wrap,stencil_type,M,N,P,m,n,p,1,s,
110:                     lx,ly,lz,&da);
111:   DAView(da,PETSC_VIEWER_DRAW_WORLD);
112:   DACreateGlobalVector(da,&global);
113:   DACreateLocalVector(da,&local);

115:   GenerateSliceScatter(da,&scatter,&vslice);

117:   /* Put the value rank+1 into all locations of vslice and transfer back to global vector */
118:   value = 1.0 + rank;
119:   VecSet(&value,vslice);
120:   VecScatterBegin(vslice,global,INSERT_VALUES,SCATTER_REVERSE,scatter);
121:   VecScatterEnd(vslice,global,INSERT_VALUES,SCATTER_REVERSE,scatter);

123:   VecView(global,PETSC_VIEWER_DRAW_WORLD);

125:   VecDestroy(local);
126:   VecDestroy(global);
127:   DADestroy(da);
128:   PetscFinalize();
129:   return 0;
130: }