Actual source code: daindex.c


  2: /*
  3:   Code for manipulating distributed regular arrays in parallel.
  4: */

  6: #include <petsc/private/dmdaimpl.h>

  8: /*
  9:    Gets the natural number for each global number on the process.

 11:    Used by DMDAGetAO() and DMDAGlobalToNatural_Create()
 12: */
 13: PetscErrorCode DMDAGetNatural_Private(DM da,PetscInt *outNlocal,IS *isnatural)
 14: {
 15:   PetscInt       Nlocal,i,j,k,*lidx,lict = 0,dim = da->dim;
 16:   DM_DA          *dd = (DM_DA*)da->data;

 18:   Nlocal = (dd->xe-dd->xs);
 19:   if (dim > 1) Nlocal *= (dd->ye-dd->ys);
 20:   if (dim > 2) Nlocal *= (dd->ze-dd->zs);

 22:   PetscMalloc1(Nlocal,&lidx);

 24:   if (dim == 1) {
 25:     for (i=dd->xs; i<dd->xe; i++) {
 26:       /*  global number in natural ordering */
 27:       lidx[lict++] = i;
 28:     }
 29:   } else if (dim == 2) {
 30:     for (j=dd->ys; j<dd->ye; j++) {
 31:       for (i=dd->xs; i<dd->xe; i++) {
 32:         /*  global number in natural ordering */
 33:         lidx[lict++] = i + j*dd->M*dd->w;
 34:       }
 35:     }
 36:   } else if (dim == 3) {
 37:     for (k=dd->zs; k<dd->ze; k++) {
 38:       for (j=dd->ys; j<dd->ye; j++) {
 39:         for (i=dd->xs; i<dd->xe; i++) {
 40:           lidx[lict++] = i + j*dd->M*dd->w + k*dd->M*dd->N*dd->w;
 41:         }
 42:       }
 43:     }
 44:   }
 45:   *outNlocal = Nlocal;
 46:   ISCreateGeneral(PetscObjectComm((PetscObject)da),Nlocal,lidx,PETSC_OWN_POINTER,isnatural);
 47:   return 0;
 48: }

 50: /*@C
 51:    DMDASetAOType - Sets the type of application ordering for a distributed array.

 53:    Collective on da

 55:    Input Parameters:
 56: +  da - the distributed array
 57: -  aotype - type of AO

 59:    Output Parameters:

 61:    Level: intermediate

 63:    Notes:
 64:    It will generate and error if an AO has already been obtained with a call to DMDAGetAO and the user sets a different AOType

 66: .seealso: DMDACreate2d(), DMDAGetAO(), DMDAGetGhostCorners(), DMDAGetCorners(), DMLocalToGlobal()
 67:           DMGlobalToLocalBegin(), DMGlobalToLocalEnd(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDAGetGlobalIndices(), DMDAGetOwnershipRanges(),
 68:           AO, AOPetscToApplication(), AOApplicationToPetsc()
 69: @*/
 70: PetscErrorCode  DMDASetAOType(DM da,AOType aotype)
 71: {
 72:   DM_DA          *dd;
 73:   PetscBool      isdmda;

 76:   PetscObjectTypeCompare((PetscObject)da,DMDA,&isdmda);
 78:   /* now we can safely dereference */
 79:   dd = (DM_DA*)da->data;
 80:   if (dd->ao) { /* check if the already computed AO has the same type as requested */
 81:     PetscBool match;
 82:     PetscObjectTypeCompare((PetscObject)dd->ao,aotype,&match);
 84:     return 0;
 85:   }
 86:   PetscFree(dd->aotype);
 87:   PetscStrallocpy(aotype,(char**)&dd->aotype);
 88:   return 0;
 89: }

 91: /*@
 92:    DMDAGetAO - Gets the application ordering context for a distributed array.

 94:    Collective on da

 96:    Input Parameter:
 97: .  da - the distributed array

 99:    Output Parameters:
100: .  ao - the application ordering context for DMDAs

102:    Level: intermediate

104:    Notes:
105:    In this case, the AO maps to the natural grid ordering that would be used
106:    for the DMDA if only 1 processor were employed (ordering most rapidly in the
107:    x-direction, then y, then z).  Multiple degrees of freedom are numbered
108:    for each node (rather than 1 component for the whole grid, then the next
109:    component, etc.)

111:    Do NOT call AODestroy() on the ao returned by this function.

113: .seealso: DMDACreate2d(), DMDASetAOType(), DMDAGetGhostCorners(), DMDAGetCorners(), DMLocalToGlobal()
114:           DMGlobalToLocalBegin(), DMGlobalToLocalEnd(), DMLocalToLocalBegin(), DMLocalToLocalEnd(),  DMDAGetOwnershipRanges(),
115:           AO, AOPetscToApplication(), AOApplicationToPetsc()
116: @*/
117: PetscErrorCode  DMDAGetAO(DM da,AO *ao)
118: {
119:   DM_DA          *dd;
120:   PetscBool      isdmda;

124:   PetscObjectTypeCompare((PetscObject)da,DMDA,&isdmda);
126:   /* now we can safely dereference */
127:   dd = (DM_DA*)da->data;

129:   /*
130:      Build the natural ordering to PETSc ordering mappings.
131:   */
132:   if (!dd->ao) {
133:     IS             ispetsc,isnatural;
134:     PetscInt       Nlocal;

136:     DMDAGetNatural_Private(da,&Nlocal,&isnatural);
137:     ISCreateStride(PetscObjectComm((PetscObject)da),Nlocal,dd->base,1,&ispetsc);
138:     AOCreate(PetscObjectComm((PetscObject)da),&dd->ao);
139:     AOSetIS(dd->ao,isnatural,ispetsc);
140:     AOSetType(dd->ao,dd->aotype);
141:     PetscLogObjectParent((PetscObject)da,(PetscObject)dd->ao);
142:     ISDestroy(&ispetsc);
143:     ISDestroy(&isnatural);
144:   }
145:   *ao = dd->ao;
146:   return 0;
147: }