Actual source code: daindex.c

  1: /*$Id: daindex.c,v 1.33 2001/06/21 21:19:09 bsmith Exp $*/
  2: 
  3: /*
  4:   Code for manipulating distributed regular arrays in parallel.
  5: */

 7:  #include src/dm/da/daimpl.h

  9: /*@C
 10:    DAGetGlobalIndices - Returns the global node number of all local nodes,
 11:    including ghost nodes.

 13:    Not Collective

 15:    Input Parameter:
 16: .  da - the distributed array

 18:    Output Parameters:
 19: +  n - the number of local elements, including ghost nodes (or PETSC_NULL)
 20: -  idx - the global indices

 22:    Level: intermediate

 24:    Note: 
 25:    For DA_STENCIL_STAR stencils the inactive corner ghost nodes are also included
 26:    in the list of local indices (even though those nodes are not updated 
 27:    during calls to DAXXXToXXX().

 29:    Essentially the same data is returned in the form of a local-to-global mapping
 30:    with the routine DAGetISLocalToGlobalMapping();

 32:    Fortran Note:
 33:    This routine is used differently from Fortran
 34: .vb
 35:         DA          da
 36:         integer     n,da_array(1)
 37:         PetscOffset i_da
 38:         integer     ierr
 39:         call DAGetGlobalIndices(da,n,da_array,i_da,ierr)

 41:    C Access first local entry in list
 42:         value = da_array(i_da + 1)
 43: .ve

 45:    See the Fortran chapter of the users manual for details

 47: .keywords: distributed array, get, global, indices, local-to-global

 49: .seealso: DACreate2d(), DAGetGhostCorners(), DAGetCorners(), DALocalToGlobal()
 50:           DAGlobalToLocal(), DALocalToLocal(), DAGetAO(), DAGetGlobalIndicesF90()
 51:           DAGetISLocalToGlobalMapping(), DACreate3d(), DACreate1d()
 52: @*/
 53: int DAGetGlobalIndices(DA da,int *n,int **idx)
 54: {
 57:   if (n)   *n   = da->Nl;
 58:   if (idx) *idx = da->idx;
 59:   return(0);
 60: }


 63: /*@C
 64:    DAGetAO - Gets the application ordering context for a distributed array.

 66:    Collective on DA

 68:    Input Parameter:
 69: .  da - the distributed array

 71:    Output Parameters:
 72: .  ao - the application ordering context for DAs

 74:    Level: intermediate

 76:    Notes:
 77:    In this case, the AO maps to the natural grid ordering that would be used
 78:    for the DA if only 1 processor were employed (ordering most rapidly in the
 79:    x-direction, then y, then z).  Multiple degrees of freedom are numbered
 80:    for each node (rather than 1 component for the whole grid, then the next
 81:    component, etc.)

 83: .keywords: distributed array, get, global, indices, local-to-global

 85: .seealso: DACreate2d(), DAGetGhostCorners(), DAGetCorners(), DALocalToGlocal()
 86:           DAGlobalToLocal(), DALocalToLocal(), DAGetGlobalIndices()
 87: @*/
 88: int DAGetAO(DA da,AO *ao)
 89: {

 93:   /* 
 94:      Build the natural ordering to PETSc ordering mappings.
 95:   */
 96:   if (!da->ao) {
 97:     IS  ispetsc,isnatural;
 98:     int ierr,i,j,k,*lidx,lict = 0,Nlocal;

100:     Nlocal = (da->xe-da->xs);
101:     if (da->dim > 1) {
102:       Nlocal *= (da->ye-da->ys);
103:     }
104:     if (da->dim > 2) {
105:       Nlocal *= (da->ze-da->zs);
106:     }

108:     ISCreateStride(da->comm,Nlocal,da->base,1,&ispetsc);
109:     PetscMalloc(Nlocal*sizeof(int),&lidx);

111:     if (da->dim == 1) {
112:        for (i=da->xs; i<da->xe; i++) {
113:          /*  global number in natural ordering */
114:          lidx[lict++] = i;
115:        }
116:     } else if (da->dim == 2) {
117:       for (j=da->ys; j<da->ye; j++) {
118:         for (i=da->xs; i<da->xe; i++) {
119:           /*  global number in natural ordering */
120:           lidx[lict++] = i + j*da->M*da->w;
121:         }
122:       }
123:     } else if (da->dim == 3) {
124:       for (k=da->zs; k<da->ze; k++) {
125:         for (j=da->ys; j<da->ye; j++) {
126:           for (i=da->xs; i<da->xe; i++) {
127:             lidx[lict++] = i + j*da->M*da->w + k*da->M*da->N*da->w;
128:           }
129:         }
130:       }
131:     }

133:     ISCreateGeneral(da->comm,Nlocal,lidx,&isnatural);
134:     PetscFree(lidx);

136:     AOCreateBasicIS(isnatural,ispetsc,&da->ao);
137:     PetscLogObjectParent(da,da->ao);
138:     ISDestroy(ispetsc);
139:     ISDestroy(isnatural);
140:   }
141:   *ao = da->ao;
142:   return(0);
143: }

145: /*MC
146:     DAGetGlobalIndicesF90 - Returns a Fortran90 pointer to the list of 
147:     global indices (global node number of all local nodes, including
148:     ghost nodes).

150:     Synopsis:
151:     DAGetGlobalIndicesF90(DA da,integer n,{integer, pointer :: idx(:)},integer ierr)

153:     Input Parameter:
154: .   da - the distributed array

156:     Output Parameters:
157: +   n - the number of local elements, including ghost nodes (or PETSC_NULL)
158: .   idx - the Fortran90 pointer to the global indices
159: -   ierr - error code

161:     Level: intermediate

163:     Notes:
164:      Not yet supported for all F90 compilers

166: .keywords: distributed array, get, global, indices, local-to-global, f90

168: .seealso: DAGetGlobalIndices()
169: M*/