Actual source code: daltol.c


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

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

  8: /*
  9:    DMLocalToLocalCreate_DA - Creates the local to local scatter

 11:    Collective on da

 13:    Input Parameter:
 14: .  da - the distributed array

 16: */
 17: PetscErrorCode  DMLocalToLocalCreate_DA(DM da)
 18: {
 19:   PetscInt       *idx,left,j,count,up,down,i,bottom,top,k,dim=da->dim;
 20:   DM_DA          *dd = (DM_DA*)da->data;


 24:   if (dd->ltol) return 0;
 25:   /*
 26:      We simply remap the values in the from part of
 27:      global to local to read from an array with the ghost values
 28:      rather then from the plain array.
 29:   */
 30:   VecScatterCopy(dd->gtol,&dd->ltol);
 31:   PetscLogObjectParent((PetscObject)da,(PetscObject)dd->ltol);
 32:   if (dim == 1) {
 33:     left = dd->xs - dd->Xs;
 34:     PetscMalloc1(dd->xe-dd->xs,&idx);
 35:     for (j=0; j<dd->xe-dd->xs; j++) idx[j] = left + j;
 36:   } else if (dim == 2) {
 37:     left  = dd->xs - dd->Xs; down  = dd->ys - dd->Ys; up    = down + dd->ye-dd->ys;
 38:     PetscMalloc1((dd->xe-dd->xs)*(up - down),&idx);
 39:     count = 0;
 40:     for (i=down; i<up; i++) {
 41:       for (j=0; j<dd->xe-dd->xs; j++) {
 42:         idx[count++] = left + i*(dd->Xe-dd->Xs) + j;
 43:       }
 44:     }
 45:   } else if (dim == 3) {
 46:     left   = dd->xs - dd->Xs;
 47:     bottom = dd->ys - dd->Ys; top = bottom + dd->ye-dd->ys;
 48:     down   = dd->zs - dd->Zs; up  = down + dd->ze-dd->zs;
 49:     count  = (dd->xe-dd->xs)*(top-bottom)*(up-down);
 50:     PetscMalloc1(count,&idx);
 51:     count  = 0;
 52:     for (i=down; i<up; i++) {
 53:       for (j=bottom; j<top; j++) {
 54:         for (k=0; k<dd->xe-dd->xs; k++) {
 55:           idx[count++] = (left+j*(dd->Xe-dd->Xs))+i*(dd->Xe-dd->Xs)*(dd->Ye-dd->Ys) + k;
 56:         }
 57:       }
 58:     }
 59:   } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_CORRUPT,"DMDA has invalid dimension %D",dim);

 61:   VecScatterRemap(dd->ltol,idx,NULL);
 62:   PetscFree(idx);
 63:   return 0;
 64: }

 66: /*
 67:    DMLocalToLocalBegin_DA - Maps from a local vector (including ghost points
 68:    that contain irrelevant values) to another local vector where the ghost
 69:    points in the second are set correctly. Must be followed by DMLocalToLocalEnd_DA().

 71:    Neighbor-wise Collective on da

 73:    Input Parameters:
 74: +  da - the distributed array context
 75: .  g - the original local vector
 76: -  mode - one of INSERT_VALUES or ADD_VALUES

 78:    Output Parameter:
 79: .  l  - the local vector with correct ghost values

 81:    Notes:
 82:    The local vectors used here need not be the same as those
 83:    obtained from DMCreateLocalVector(), BUT they
 84:    must have the same parallel data layout; they could, for example, be
 85:    obtained with VecDuplicate() from the DMDA originating vectors.

 87: */
 88: PetscErrorCode  DMLocalToLocalBegin_DA(DM da,Vec g,InsertMode mode,Vec l)
 89: {
 90:   DM_DA          *dd = (DM_DA*)da->data;

 93:   if (!dd->ltol) {
 94:     DMLocalToLocalCreate_DA(da);
 95:   }
 96:   VecScatterBegin(dd->ltol,g,l,mode,SCATTER_FORWARD);
 97:   return 0;
 98: }

100: /*
101:    DMLocalToLocalEnd_DA - Maps from a local vector (including ghost points
102:    that contain irrelevant values) to another local vector where the ghost
103:    points in the second are set correctly.  Must be preceded by
104:    DMLocalToLocalBegin_DA().

106:    Neighbor-wise Collective on da

108:    Input Parameters:
109: +  da - the distributed array context
110: .  g - the original local vector
111: -  mode - one of INSERT_VALUES or ADD_VALUES

113:    Output Parameter:
114: .  l  - the local vector with correct ghost values

116:    Note:
117:    The local vectors used here need not be the same as those
118:    obtained from DMCreateLocalVector(), BUT they
119:    must have the same parallel data layout; they could, for example, be
120:    obtained with VecDuplicate() from the DMDA originating vectors.

122: */
123: PetscErrorCode  DMLocalToLocalEnd_DA(DM da,Vec g,InsertMode mode,Vec l)
124: {
125:   DM_DA          *dd = (DM_DA*)da->data;

129:   VecScatterEnd(dd->ltol,g,l,mode,SCATTER_FORWARD);
130:   return 0;
131: }