Actual source code: daltol.c

petsc-3.3-p7 2013-05-11
  2: /*
  3:   Code for manipulating distributed regular arrays in parallel.
  4: */

  6: #include <petsc-private/daimpl.h>    /*I   "petscdmda.h"   I*/

 10: /*
 11:    DMDALocalToLocalCreate - Creates the local to local scatter

 13:    Collective on DMDA

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

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


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

 67:   VecScatterRemap(dd->ltol,idx,PETSC_NULL);
 68:   PetscFree(idx);
 69:   return(0);
 70: }

 74: /*@
 75:    DMDALocalToLocalBegin - Maps from a local vector (including ghost points
 76:    that contain irrelevant values) to another local vector where the ghost
 77:    points in the second are set correctly. Must be followed by DMDALocalToLocalEnd().

 79:    Neighbor-wise Collective on DMDA and Vec

 81:    Input Parameters:
 82: +  da - the distributed array context
 83: .  g - the original local vector
 84: -  mode - one of INSERT_VALUES or ADD_VALUES

 86:    Output Parameter:
 87: .  l  - the local vector with correct ghost values

 89:    Level: intermediate

 91:    Notes:
 92:    The local vectors used here need not be the same as those
 93:    obtained from DMCreateLocalVector(), BUT they
 94:    must have the same parallel data layout; they could, for example, be 
 95:    obtained with VecDuplicate() from the DMDA originating vectors.

 97: .keywords: distributed array, local-to-local, begin

 99: .seealso: DMDALocalToLocalEnd(), DMLocalToGlobalBegin()
100: @*/
101: PetscErrorCode  DMDALocalToLocalBegin(DM da,Vec g,InsertMode mode,Vec l)
102: {
104:   DM_DA          *dd = (DM_DA*)da->data;

108:   if (!dd->ltol) {
109:     DMDALocalToLocalCreate(da);
110:   }
111:   VecScatterBegin(dd->ltol,g,l,mode,SCATTER_FORWARD);
112:   return(0);
113: }

117: /*@
118:    DMDALocalToLocalEnd - Maps from a local vector (including ghost points
119:    that contain irrelevant values) to another local vector where the ghost
120:    points in the second are set correctly.  Must be preceeded by 
121:    DMDALocalToLocalBegin().

123:    Neighbor-wise Collective on DMDA and Vec

125:    Input Parameters:
126: +  da - the distributed array context
127: .  g - the original local vector
128: -  mode - one of INSERT_VALUES or ADD_VALUES

130:    Output Parameter:
131: .  l  - the local vector with correct ghost values

133:    Level: intermediate

135:    Note:
136:    The local vectors used here need not be the same as those
137:    obtained from DMCreateLocalVector(), BUT they
138:    must have the same parallel data layout; they could, for example, be 
139:    obtained with VecDuplicate() from the DMDA originating vectors.

141: .keywords: distributed array, local-to-local, end

143: .seealso: DMDALocalToLocalBegin(), DMLocalToGlobalBegin()
144: @*/
145: PetscErrorCode  DMDALocalToLocalEnd(DM da,Vec g,InsertMode mode,Vec l)
146: {
148:   DM_DA          *dd = (DM_DA*)da->data;

154:   VecScatterEnd(dd->ltol,g,l,mode,SCATTER_FORWARD);
155:   return(0);
156: }