Actual source code: da3.c

  1: /*
  2:    Code for manipulating distributed regular 3d arrays in parallel.
  3:    File created by Peter Mell  7/14/95
  4:  */

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

 10: PetscErrorCode DAView_3d(DA da,PetscViewer viewer)
 11: {
 13:   PetscMPIInt    rank;
 14:   PetscTruth     iascii,isdraw;

 17:   MPI_Comm_rank(da->comm,&rank);

 19:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
 20:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
 21:   if (iascii) {
 22:     PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %D N %D P %D m %D n %D p %D w %D s %D\n",
 23:                rank,da->M,da->N,da->P,da->m,da->n,da->p,da->w,da->s);
 24:     PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D, Y range of indices: %D %D, Z range of indices: %D %D\n",
 25:                da->xs,da->xe,da->ys,da->ye,da->zs,da->ze);
 26: #if !defined(PETSC_USE_COMPLEX)
 27:     if (da->coordinates) {
 28:       PetscInt  last;
 29:       PetscReal *coors;
 30:       VecGetArray(da->coordinates,&coors);
 31:       VecGetLocalSize(da->coordinates,&last);
 32:       last = last - 3;
 33:       PetscViewerASCIISynchronizedPrintf(viewer,"Lower left corner %g %g %g : Upper right %g %g %g\n",
 34:                coors[0],coors[1],coors[2],coors[last],coors[last+1],coors[last+2]);
 35:       VecRestoreArray(da->coordinates,&coors);
 36:     }
 37: #endif
 38:     PetscViewerFlush(viewer);
 39:   } else if (isdraw) {
 40:     PetscDraw       draw;
 41:     PetscReal     ymin = -1.0,ymax = (PetscReal)da->N;
 42:     PetscReal     xmin = -1.0,xmax = (PetscReal)((da->M+2)*da->P),x,y,ycoord,xcoord;
 43:     PetscInt        k,plane,base,*idx;
 44:     char       node[10];
 45:     PetscTruth isnull;

 47:     PetscViewerDrawGetDraw(viewer,0,&draw);
 48:     PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
 49:     PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);
 50:     PetscDrawSynchronizedClear(draw);

 52:     /* first processor draw all node lines */
 53:     if (!rank) {
 54:       for (k=0; k<da->P; k++) {
 55:         ymin = 0.0; ymax = (PetscReal)(da->N - 1);
 56:         for (xmin=(PetscReal)(k*(da->M+1)); xmin<(PetscReal)(da->M+(k*(da->M+1))); xmin++) {
 57:           PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);
 58:         }
 59: 
 60:         xmin = (PetscReal)(k*(da->M+1)); xmax = xmin + (PetscReal)(da->M - 1);
 61:         for (ymin=0; ymin<(PetscReal)da->N; ymin++) {
 62:           PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);
 63:         }
 64:       }
 65:     }
 66:     PetscDrawSynchronizedFlush(draw);
 67:     PetscDrawPause(draw);

 69:     for (k=0; k<da->P; k++) {  /*Go through and draw for each plane*/
 70:       if ((k >= da->zs) && (k < da->ze)) {
 71:         /* draw my box */
 72:         ymin = da->ys;
 73:         ymax = da->ye - 1;
 74:         xmin = da->xs/da->w    + (da->M+1)*k;
 75:         xmax =(da->xe-1)/da->w + (da->M+1)*k;

 77:         PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);
 78:         PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);
 79:         PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);
 80:         PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);

 82:         xmin = da->xs/da->w;
 83:         xmax =(da->xe-1)/da->w;

 85:         /* put in numbers*/
 86:         base = (da->base+(da->xe-da->xs)*(da->ye-da->ys)*(k-da->zs))/da->w;

 88:         /* Identify which processor owns the box */
 89:         sprintf(node,"%d",rank);
 90:         PetscDrawString(draw,xmin+(da->M+1)*k+.2,ymin+.3,PETSC_DRAW_RED,node);

 92:         for (y=ymin; y<=ymax; y++) {
 93:           for (x=xmin+(da->M+1)*k; x<=xmax+(da->M+1)*k; x++) {
 94:             sprintf(node,"%d",(int)base++);
 95:             PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);
 96:           }
 97:         }
 98: 
 99:       }
100:     }
101:     PetscDrawSynchronizedFlush(draw);
102:     PetscDrawPause(draw);

104:     for (k=0-da->s; k<da->P+da->s; k++) {
105:       /* Go through and draw for each plane */
106:       if ((k >= da->Zs) && (k < da->Ze)) {
107: 
108:         /* overlay ghost numbers, useful for error checking */
109:         base = (da->Xe-da->Xs)*(da->Ye-da->Ys)*(k-da->Zs); idx = da->idx;
110:         plane=k;
111:         /* Keep z wrap around points on the dradrawg */
112:         if (k<0)    { plane=da->P+k; }
113:         if (k>=da->P) { plane=k-da->P; }
114:         ymin = da->Ys; ymax = da->Ye;
115:         xmin = (da->M+1)*plane*da->w;
116:         xmax = (da->M+1)*plane*da->w+da->M*da->w;
117:         for (y=ymin; y<ymax; y++) {
118:           for (x=xmin+da->Xs; x<xmin+da->Xe; x+=da->w) {
119:             sprintf(node,"%d",(int)(idx[base]/da->w));
120:             ycoord = y;
121:             /*Keep y wrap around points on drawing */
122:             if (y<0)      { ycoord = da->N+y; }

124:             if (y>=da->N) { ycoord = y-da->N; }
125:             xcoord = x;   /* Keep x wrap points on drawing */

127:             if (x<xmin)  { xcoord = xmax - (xmin-x); }
128:             if (x>=xmax) { xcoord = xmin + (x-xmax); }
129:             PetscDrawString(draw,xcoord/da->w,ycoord,PETSC_DRAW_BLUE,node);
130:             base+=da->w;
131:           }
132:         }
133:       }
134:     }
135:     PetscDrawSynchronizedFlush(draw);
136:     PetscDrawPause(draw);
137:   } else {
138:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for DA 3d",((PetscObject)viewer)->type_name);
139:   }
140:   return(0);
141: }

143: EXTERN PetscErrorCode DAPublish_Petsc(PetscObject);

147: /*@C
148:    DACreate3d - Creates an object that will manage the communication of three-dimensional 
149:    regular array data that is distributed across some processors.

151:    Collective on MPI_Comm

153:    Input Parameters:
154: +  comm - MPI communicator
155: .  wrap - type of periodicity the array should have, if any.  Use one
156:           of DA_NONPERIODIC, DA_XPERIODIC, DA_YPERIODIC, DA_XYPERIODIC, DA_XYZPERIODIC, DA_XZPERIODIC, or DA_YZPERIODIC.
157: .  stencil_type - Type of stencil (DA_STENCIL_STAR or DA_STENCIL_BOX)
158: .  M,N,P - global dimension in each direction of the array (use -M, -N, and or -P to indicate that it may be set to a different value 
159:             from the command line with -da_grid_x <M> -da_grid_y <N> -da_grid_z <P>)
160: .  m,n,p - corresponding number of processors in each dimension 
161:            (or PETSC_DECIDE to have calculated)
162: .  dof - number of degrees of freedom per node
163: .  lx, ly, lz - arrays containing the number of nodes in each cell along
164:           the x, y, and z coordinates, or PETSC_NULL. If non-null, these
165:           must be of length as m,n,p and the corresponding
166:           m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of
167:           the ly[] must n, sum of the lz[] must be P
168: -  s - stencil width

170:    Output Parameter:
171: .  inra - the resulting distributed array object

173:    Options Database Key:
174: +  -da_view - Calls DAView() at the conclusion of DACreate3d()
175: .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
176: .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
177: .  -da_grid_z <nz> - number of grid points in z direction, if P < 0
178: .  -da_refine_x - refinement ratio in x direction
179: .  -da_refine_y - refinement ratio in y direction
180: -  -da_refine_y - refinement ratio in z direction

182:    Level: beginner

184:    Notes:
185:    The stencil type DA_STENCIL_STAR with width 1 corresponds to the 
186:    standard 7-pt stencil, while DA_STENCIL_BOX with width 1 denotes
187:    the standard 27-pt stencil.

189:    The array data itself is NOT stored in the DA, it is stored in Vec objects;
190:    The appropriate vector objects can be obtained with calls to DACreateGlobalVector()
191:    and DACreateLocalVector() and calls to VecDuplicate() if more are needed.

193: .keywords: distributed array, create, three-dimensional

195: .seealso: DADestroy(), DAView(), DACreate1d(), DACreate2d(), DAGlobalToLocalBegin(), DAGetRefinementFactor(),
196:           DAGlobalToLocalEnd(), DALocalToGlobal(), DALocalToLocalBegin(), DALocalToLocalEnd(), DASetRefinementFactor(),
197:           DAGetInfo(), DACreateGlobalVector(), DACreateLocalVector(), DACreateNaturalVector(), DALoad(), DAView()

199: @*/
200: PetscErrorCode DACreate3d(MPI_Comm comm,DAPeriodicType wrap,DAStencilType stencil_type,PetscInt M,
201:                PetscInt N,PetscInt P,PetscInt m,PetscInt n,PetscInt p,PetscInt dof,PetscInt s,PetscInt *lx,PetscInt *ly,PetscInt *lz,DA *inra)
202: {
204:   PetscMPIInt    rank,size;
205:   PetscInt       xs,xe,ys,ye,zs,ze,x,y,z,Xs,Xe,Ys,Ye,Zs,Ze,start,end,pm;
206:   PetscInt       left,up,down,bottom,top,i,j,k,*idx,nn,*flx = 0,*fly = 0,*flz = 0;
207:   PetscInt       n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14;
208:   PetscInt       n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26;
209:   PetscInt       *bases,*ldims,x_t,y_t,z_t,s_t,base,count,s_x,s_y,s_z;
210:   PetscInt       tM = M,tN = N,tP = P;
211:   PetscInt       sn0 = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0;
212:   PetscInt       sn8 = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0;
213:   PetscInt       sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0,refine_x = 2, refine_y = 2, refine_z = 2;
214:   PetscTruth     flg1,flg2;
215:   DA             da;
216:   Vec            local,global;
217:   VecScatter     ltog,gtol;
218:   IS             to,from;

222:   *inra = 0;
223: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
224:   DMInitializePackage(PETSC_NULL);
225: #endif

227:   if (dof < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
228:   if (s < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);

230:   PetscOptionsBegin(comm,PETSC_NULL,"3d DA Options","DA");
231:     if (M < 0){
232:       tM   = -M;
233:       PetscOptionsInt("-da_grid_x","Number of grid points in x direction","DACreate3d",tM,&tM,PETSC_NULL);
234:     }
235:     if (N < 0){
236:       tN   = -N;
237:       PetscOptionsInt("-da_grid_y","Number of grid points in y direction","DACreate3d",tN,&tN,PETSC_NULL);
238:     }
239:     if (P < 0){
240:       tP   = -P;
241:       PetscOptionsInt("-da_grid_z","Number of grid points in z direction","DACreate3d",tP,&tP,PETSC_NULL);
242:     }
243:     PetscOptionsInt("-da_processors_x","Number of processors in x direction","DACreate3d",m,&m,PETSC_NULL);
244:     PetscOptionsInt("-da_processors_y","Number of processors in y direction","DACreate3d",n,&n,PETSC_NULL);
245:     PetscOptionsInt("-da_processors_z","Number of processors in z direction","DACreate3d",p,&p,PETSC_NULL);
246:     PetscOptionsInt("-da_refine_x","Refinement ratio in x direction","DASetRefinementFactor",refine_x,&refine_x,PETSC_NULL);
247:     PetscOptionsInt("-da_refine_y","Refinement ratio in y direction","DASetRefinementFactor",refine_y,&refine_y,PETSC_NULL);
248:     PetscOptionsInt("-da_refine_z","Refinement ratio in z direction","DASetRefinementFactor",refine_z,&refine_z,PETSC_NULL);
249:   PetscOptionsEnd();
250:   M = tM; N = tN; P = tP;

252:   PetscHeaderCreate(da,_p_DA,struct _DAOps,DA_COOKIE,0,"DA",comm,DADestroy,DAView);
253:   da->bops->publish           = DAPublish_Petsc;
254:   da->ops->createglobalvector = DACreateGlobalVector;
255:   da->ops->getinterpolation   = DAGetInterpolation;
256:   da->ops->getcoloring        = DAGetColoring;
257:   da->ops->getmatrix          = DAGetMatrix;
258:   da->ops->refine             = DARefine;

260:   PetscLogObjectCreate(da);
261:   PetscLogObjectMemory(da,sizeof(struct _p_DA));
262:   da->dim        = 3;
263:   da->interptype = DA_Q1;
264:   da->refine_x   = refine_x;
265:   da->refine_y   = refine_y;
266:   da->refine_z   = refine_z;
267:   PetscMalloc(dof*sizeof(char*),&da->fieldname);
268:   PetscMemzero(da->fieldname,dof*sizeof(char*));

270:   MPI_Comm_size(comm,&size);
271:   MPI_Comm_rank(comm,&rank);

273:   if (m != PETSC_DECIDE) {
274:     if (m < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);}
275:     else if (m > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);}
276:   }
277:   if (n != PETSC_DECIDE) {
278:     if (n < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);}
279:     else if (n > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);}
280:   }
281:   if (p != PETSC_DECIDE) {
282:     if (p < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %D",p);}
283:     else if (p > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %D %d",p,size);}
284:   }

286:   /* Partition the array among the processors */
287:   if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) {
288:     m = size/(n*p);
289:   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
290:     n = size/(m*p);
291:   } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
292:     p = size/(m*n);
293:   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
294:     /* try for squarish distribution */
295:     m = (int)(0.5 + sqrt(((PetscReal)M)*((PetscReal)size)/((PetscReal)N*p)));
296:     if (!m) m = 1;
297:     while (m > 0) {
298:       n = size/(m*p);
299:       if (m*n*p == size) break;
300:       m--;
301:     }
302:     if (!m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %D",p);
303:     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
304:   } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
305:     /* try for squarish distribution */
306:     m = (int)(0.5 + sqrt(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n)));
307:     if (!m) m = 1;
308:     while (m > 0) {
309:       p = size/(m*n);
310:       if (m*n*p == size) break;
311:       m--;
312:     }
313:     if (!m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %D",n);
314:     if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;}
315:   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
316:     /* try for squarish distribution */
317:     n = (int)(0.5 + sqrt(((PetscReal)N)*((PetscReal)size)/((PetscReal)P*m)));
318:     if (!n) n = 1;
319:     while (n > 0) {
320:       p = size/(m*n);
321:       if (m*n*p == size) break;
322:       n--;
323:     }
324:     if (!n) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %D",n);
325:     if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;}
326:   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
327:     /* try for squarish distribution */
328:     n = (PetscInt)(0.5 + pow(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),1./3.));
329:     if (!n) n = 1;
330:     while (n > 0) {
331:       pm = size/n;
332:       if (n*pm == size) break;
333:       n--;
334:     }
335:     if (!n) n = 1;
336:     m = (PetscInt)(0.5 + sqrt(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n)));
337:     if (!m) m = 1;
338:     while (m > 0) {
339:       p = size/(m*n);
340:       if (m*n*p == size) break;
341:       m--;
342:     }
343:     if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;}
344:   } else if (m*n*p != size) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");

346:   if (m*n*p != size) SETERRQ(PETSC_ERR_PLIB,"Could not find good partition");
347:   if (M < m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
348:   if (N < n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
349:   if (P < p) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %D %D",P,p);

351:   PetscOptionsHasName(PETSC_NULL,"-da_partition_nodes_at_end",&flg2);
352:   /* 
353:      Determine locally owned region 
354:      [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes 
355:   */
356:   if (lx) { /* user decided distribution */
357:     x  = lx[rank % m];
358:     xs = 0;
359:     for (i=0; i<(rank%m); i++) { xs += lx[i];}
360:     if (m > 1 && x < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column width is too thin for stencil! %D %D",x,s);
361:   } else if (flg2) {
362:     SETERRQ(PETSC_ERR_SUP,"-da_partition_nodes_at_end not supported");
363:   } else { /* Normal PETSc distribution */
364:     x = M/m + ((M % m) > (rank % m));
365:     if (m > 1 && x < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column width is too thin for stencil! %D %D",x,s);
366:     if ((M % m) > (rank % m)) { xs = (rank % m)*x; }
367:     else                      { xs = (M % m)*(x+1) + ((rank % m)-(M % m))*x; }
368:     PetscMalloc(m*sizeof(PetscInt),&lx);
369:     flx = lx;
370:     for (i=0; i<m; i++) {
371:       lx[i] = M/m + ((M % m) > (i % m));
372:     }
373:   }
374:   if (ly) { /* user decided distribution */
375:     y  = ly[(rank % (m*n))/m];
376:     if (n > 1 && y < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row width is too thin for stencil! %D %D",y,s);
377:     ys = 0;
378:     for (i=0; i<(rank % (m*n))/m; i++) { ys += ly[i];}
379:   } else if (flg2) {
380:     SETERRQ(PETSC_ERR_SUP,"-da_partition_nodes_at_end not supported");
381:   } else { /* Normal PETSc distribution */
382:     y = N/n + ((N % n) > ((rank % (m*n)) /m));
383:     if (n > 1 && y < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row width is too thin for stencil! %D %D",y,s);
384:     if ((N % n) > ((rank % (m*n)) /m)) {ys = ((rank % (m*n))/m)*y;}
385:     else                               {ys = (N % n)*(y+1) + (((rank % (m*n))/m)-(N % n))*y;}
386:     PetscMalloc(n*sizeof(PetscInt),&ly);
387:     fly = ly;
388:     for (i=0; i<n; i++) {
389:       ly[i] = N/n + ((N % n) > (i % n));
390:     }
391:   }
392:   if (lz) { /* user decided distribution */
393:     z  = lz[rank/(m*n)];
394:     if (p > 1 && z < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Plane width is too thin for stencil! %D %D",z,s);
395:     zs = 0;
396:     for (i=0; i<(rank/(m*n)); i++) { zs += lz[i];}
397:   } else if (flg2) {
398:     SETERRQ(PETSC_ERR_SUP,"-da_partition_nodes_at_end not supported");
399:   } else { /* Normal PETSc distribution */
400:     z = P/p + ((P % p) > (rank / (m*n)));
401:     if (p > 1 && z < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Plane width is too thin for stencil! %D %D",z,s);
402:     if ((P % p) > (rank / (m*n))) {zs = (rank/(m*n))*z;}
403:     else                          {zs = (P % p)*(z+1) + ((rank/(m*n))-(P % p))*z;}
404:     PetscMalloc(p*sizeof(PetscInt),&lz);
405:     flz = lz;
406:     for (i=0; i<p; i++) {
407:       lz[i] = P/p + ((P % p) > (i % p));
408:     }
409:   }
410:   ye = ys + y;
411:   xe = xs + x;
412:   ze = zs + z;

414:   /* determine ghost region */
415:   /* Assume No Periodicity */
416:   if (xs-s > 0) Xs = xs - s; else Xs = 0;
417:   if (ys-s > 0) Ys = ys - s; else Ys = 0;
418:   if (zs-s > 0) Zs = zs - s; else Zs = 0;
419:   if (xe+s <= M) Xe = xe + s; else Xe = M;
420:   if (ye+s <= N) Ye = ye + s; else Ye = N;
421:   if (ze+s <= P) Ze = ze + s; else Ze = P;

423:   /* X Periodic */
424:   if (DAXPeriodic(wrap)){
425:     Xs = xs - s;
426:     Xe = xe + s;
427:   }

429:   /* Y Periodic */
430:   if (DAYPeriodic(wrap)){
431:     Ys = ys - s;
432:     Ye = ye + s;
433:   }

435:   /* Z Periodic */
436:   if (DAZPeriodic(wrap)){
437:     Zs = zs - s;
438:     Ze = ze + s;
439:   }

441:   /* Resize all X parameters to reflect w */
442:   x   *= dof;
443:   xs  *= dof;
444:   xe  *= dof;
445:   Xs  *= dof;
446:   Xe  *= dof;
447:   s_x  = s*dof;
448:   s_y  = s;
449:   s_z  = s;

451:   /* determine starting point of each processor */
452:   nn       = x*y*z;
453:   PetscMalloc((2*size+1)*sizeof(PetscInt),&bases);
454:   ldims    = (PetscInt*)(bases+size+1);
455:   MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);
456:   bases[0] = 0;
457:   for (i=1; i<=size; i++) {
458:     bases[i] = ldims[i-1];
459:   }
460:   for (i=1; i<=size; i++) {
461:     bases[i] += bases[i-1];
462:   }

464:   /* allocate the base parallel and sequential vectors */
465:   da->Nlocal = x*y*z;
466:   VecCreateMPIWithArray(comm,da->Nlocal,PETSC_DECIDE,0,&global);
467:   VecSetBlockSize(global,dof);
468:   da->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs);
469:   VecCreateSeqWithArray(MPI_COMM_SELF,da->nlocal,0,&local);
470:   VecSetBlockSize(local,dof);

472:   /* generate appropriate vector scatters */
473:   /* local to global inserts non-ghost point region into global */
474:   VecGetOwnershipRange(global,&start,&end);
475:   ISCreateStride(comm,x*y*z,start,1,&to);

477:   left   = xs - Xs;
478:   bottom = ys - Ys; top = bottom + y;
479:   down   = zs - Zs; up  = down + z;
480:   count  = x*(top-bottom)*(up-down);
481:   PetscMalloc(count*sizeof(PetscInt),&idx);
482:   count  = 0;
483:   for (i=down; i<up; i++) {
484:     for (j=bottom; j<top; j++) {
485:       for (k=0; k<x; k++) {
486:         idx[count++] = (left+j*(Xe-Xs))+i*(Xe-Xs)*(Ye-Ys) + k;
487:       }
488:     }
489:   }
490:   ISCreateGeneral(comm,count,idx,&from);
491:   PetscFree(idx);

493:   VecScatterCreate(local,from,global,to,&ltog);
494:   PetscLogObjectParent(da,to);
495:   PetscLogObjectParent(da,from);
496:   PetscLogObjectParent(da,ltog);
497:   ISDestroy(from);
498:   ISDestroy(to);

500:   /* global to local must include ghost points */
501:   if (stencil_type == DA_STENCIL_BOX) {
502:     ISCreateStride(comm,(Xe-Xs)*(Ye-Ys)*(Ze-Zs),0,1,&to);
503:   } else {
504:     /* This is way ugly! We need to list the funny cross type region */
505:     /* the bottom chunck */
506:     left   = xs - Xs;
507:     bottom = ys - Ys; top = bottom + y;
508:     down   = zs - Zs;   up  = down + z;
509:     count  = down*(top-bottom)*x +
510:              (up-down)*(bottom*x  + (top-bottom)*(Xe-Xs) + (Ye-Ys-top)*x) +
511:              (Ze-Zs-up)*(top-bottom)*x;
512:     PetscMalloc(count*sizeof(PetscInt),&idx);
513:     count  = 0;
514:     for (i=0; i<down; i++) {
515:       for (j=bottom; j<top; j++) {
516:         for (k=0; k<x; k++) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
517:       }
518:     }
519:     /* the middle piece */
520:     for (i=down; i<up; i++) {
521:       /* front */
522:       for (j=0; j<bottom; j++) {
523:         for (k=0; k<x; k++) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
524:       }
525:       /* middle */
526:       for (j=bottom; j<top; j++) {
527:         for (k=0; k<Xe-Xs; k++) idx[count++] = j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
528:       }
529:       /* back */
530:       for (j=top; j<Ye-Ys; j++) {
531:         for (k=0; k<x; k++) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
532:       }
533:     }
534:     /* the top piece */
535:     for (i=up; i<Ze-Zs; i++) {
536:       for (j=bottom; j<top; j++) {
537:         for (k=0; k<x; k++) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
538:       }
539:     }
540:     ISCreateGeneral(comm,count,idx,&to);
541:     PetscFree(idx);
542:   }

544:   /* determine who lies on each side of use stored in    n24 n25 n26
545:                                                          n21 n22 n23
546:                                                          n18 n19 n20

548:                                                          n15 n16 n17
549:                                                          n12     n14
550:                                                          n9  n10 n11

552:                                                          n6  n7  n8
553:                                                          n3  n4  n5
554:                                                          n0  n1  n2
555:   */
556: 
557:   /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
558: 
559:   /* Assume Nodes are Internal to the Cube */
560: 
561:   n0  = rank - m*n - m - 1;
562:   n1  = rank - m*n - m;
563:   n2  = rank - m*n - m + 1;
564:   n3  = rank - m*n -1;
565:   n4  = rank - m*n;
566:   n5  = rank - m*n + 1;
567:   n6  = rank - m*n + m - 1;
568:   n7  = rank - m*n + m;
569:   n8  = rank - m*n + m + 1;

571:   n9  = rank - m - 1;
572:   n10 = rank - m;
573:   n11 = rank - m + 1;
574:   n12 = rank - 1;
575:   n14 = rank + 1;
576:   n15 = rank + m - 1;
577:   n16 = rank + m;
578:   n17 = rank + m + 1;

580:   n18 = rank + m*n - m - 1;
581:   n19 = rank + m*n - m;
582:   n20 = rank + m*n - m + 1;
583:   n21 = rank + m*n - 1;
584:   n22 = rank + m*n;
585:   n23 = rank + m*n + 1;
586:   n24 = rank + m*n + m - 1;
587:   n25 = rank + m*n + m;
588:   n26 = rank + m*n + m + 1;

590:   /* Assume Pieces are on Faces of Cube */

592:   if (xs == 0) { /* First assume not corner or edge */
593:     n0  = rank       -1 - (m*n);
594:     n3  = rank + m   -1 - (m*n);
595:     n6  = rank + 2*m -1 - (m*n);
596:     n9  = rank       -1;
597:     n12 = rank + m   -1;
598:     n15 = rank + 2*m -1;
599:     n18 = rank       -1 + (m*n);
600:     n21 = rank + m   -1 + (m*n);
601:     n24 = rank + 2*m -1 + (m*n);
602:    }

604:   if (xe == M*dof) { /* First assume not corner or edge */
605:     n2  = rank -2*m +1 - (m*n);
606:     n5  = rank - m  +1 - (m*n);
607:     n8  = rank      +1 - (m*n);
608:     n11 = rank -2*m +1;
609:     n14 = rank - m  +1;
610:     n17 = rank      +1;
611:     n20 = rank -2*m +1 + (m*n);
612:     n23 = rank - m  +1 + (m*n);
613:     n26 = rank      +1 + (m*n);
614:   }

616:   if (ys==0) { /* First assume not corner or edge */
617:     n0  = rank + m * (n-1) -1 - (m*n);
618:     n1  = rank + m * (n-1)    - (m*n);
619:     n2  = rank + m * (n-1) +1 - (m*n);
620:     n9  = rank + m * (n-1) -1;
621:     n10 = rank + m * (n-1);
622:     n11 = rank + m * (n-1) +1;
623:     n18 = rank + m * (n-1) -1 + (m*n);
624:     n19 = rank + m * (n-1)    + (m*n);
625:     n20 = rank + m * (n-1) +1 + (m*n);
626:   }

628:   if (ye == N) { /* First assume not corner or edge */
629:     n6  = rank - m * (n-1) -1 - (m*n);
630:     n7  = rank - m * (n-1)    - (m*n);
631:     n8  = rank - m * (n-1) +1 - (m*n);
632:     n15 = rank - m * (n-1) -1;
633:     n16 = rank - m * (n-1);
634:     n17 = rank - m * (n-1) +1;
635:     n24 = rank - m * (n-1) -1 + (m*n);
636:     n25 = rank - m * (n-1)    + (m*n);
637:     n26 = rank - m * (n-1) +1 + (m*n);
638:   }
639: 
640:   if (zs == 0) { /* First assume not corner or edge */
641:     n0 = size - (m*n) + rank - m - 1;
642:     n1 = size - (m*n) + rank - m;
643:     n2 = size - (m*n) + rank - m + 1;
644:     n3 = size - (m*n) + rank - 1;
645:     n4 = size - (m*n) + rank;
646:     n5 = size - (m*n) + rank + 1;
647:     n6 = size - (m*n) + rank + m - 1;
648:     n7 = size - (m*n) + rank + m ;
649:     n8 = size - (m*n) + rank + m + 1;
650:   }

652:   if (ze == P) { /* First assume not corner or edge */
653:     n18 = (m*n) - (size-rank) - m - 1;
654:     n19 = (m*n) - (size-rank) - m;
655:     n20 = (m*n) - (size-rank) - m + 1;
656:     n21 = (m*n) - (size-rank) - 1;
657:     n22 = (m*n) - (size-rank);
658:     n23 = (m*n) - (size-rank) + 1;
659:     n24 = (m*n) - (size-rank) + m - 1;
660:     n25 = (m*n) - (size-rank) + m;
661:     n26 = (m*n) - (size-rank) + m + 1;
662:   }

664:   if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */
665:     n0 = size - m*n + rank + m-1 - m;
666:     n3 = size - m*n + rank + m-1;
667:     n6 = size - m*n + rank + m-1 + m;
668:   }
669: 
670:   if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */
671:     n18 = m*n - (size - rank) + m-1 - m;
672:     n21 = m*n - (size - rank) + m-1;
673:     n24 = m*n - (size - rank) + m-1 + m;
674:   }

676:   if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */
677:     n0  = rank + m*n -1 - m*n;
678:     n9  = rank + m*n -1;
679:     n18 = rank + m*n -1 + m*n;
680:   }

682:   if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */
683:     n6  = rank - m*(n-1) + m-1 - m*n;
684:     n15 = rank - m*(n-1) + m-1;
685:     n24 = rank - m*(n-1) + m-1 + m*n;
686:   }

688:   if ((xe==M*dof) && (zs==0)) { /* Assume an edge, not corner */
689:     n2 = size - (m*n-rank) - (m-1) - m;
690:     n5 = size - (m*n-rank) - (m-1);
691:     n8 = size - (m*n-rank) - (m-1) + m;
692:   }

694:   if ((xe==M*dof) && (ze==P)) { /* Assume an edge, not corner */
695:     n20 = m*n - (size - rank) - (m-1) - m;
696:     n23 = m*n - (size - rank) - (m-1);
697:     n26 = m*n - (size - rank) - (m-1) + m;
698:   }

700:   if ((xe==M*dof) && (ys==0)) { /* Assume an edge, not corner */
701:     n2  = rank + m*(n-1) - (m-1) - m*n;
702:     n11 = rank + m*(n-1) - (m-1);
703:     n20 = rank + m*(n-1) - (m-1) + m*n;
704:   }

706:   if ((xe==M*dof) && (ye==N)) { /* Assume an edge, not corner */
707:     n8  = rank - m*n +1 - m*n;
708:     n17 = rank - m*n +1;
709:     n26 = rank - m*n +1 + m*n;
710:   }

712:   if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */
713:     n0 = size - m + rank -1;
714:     n1 = size - m + rank;
715:     n2 = size - m + rank +1;
716:   }

718:   if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */
719:     n18 = m*n - (size - rank) + m*(n-1) -1;
720:     n19 = m*n - (size - rank) + m*(n-1);
721:     n20 = m*n - (size - rank) + m*(n-1) +1;
722:   }

724:   if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */
725:     n6 = size - (m*n-rank) - m * (n-1) -1;
726:     n7 = size - (m*n-rank) - m * (n-1);
727:     n8 = size - (m*n-rank) - m * (n-1) +1;
728:   }

730:   if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */
731:     n24 = rank - (size-m) -1;
732:     n25 = rank - (size-m);
733:     n26 = rank - (size-m) +1;
734:   }

736:   /* Check for Corners */
737:   if ((xs==0)   && (ys==0) && (zs==0)) { n0  = size -1;}
738:   if ((xs==0)   && (ys==0) && (ze==P)) { n18 = m*n-1;}
739:   if ((xs==0)   && (ye==N) && (zs==0)) { n6  = (size-1)-m*(n-1);}
740:   if ((xs==0)   && (ye==N) && (ze==P)) { n24 = m-1;}
741:   if ((xe==M*dof) && (ys==0) && (zs==0)) { n2  = size-m;}
742:   if ((xe==M*dof) && (ys==0) && (ze==P)) { n20 = m*n-m;}
743:   if ((xe==M*dof) && (ye==N) && (zs==0)) { n8  = size-m*n;}
744:   if ((xe==M*dof) && (ye==N) && (ze==P)) { n26 = 0;}

746:   /* Check for when not X,Y, and Z Periodic */

748:   /* If not X periodic */
749:   if ((wrap != DA_XPERIODIC)  && (wrap != DA_XYPERIODIC) &&
750:      (wrap != DA_XZPERIODIC) && (wrap != DA_XYZPERIODIC)) {
751:     if (xs==0)   {n0  = n3  = n6  = n9  = n12 = n15 = n18 = n21 = n24 = -2;}
752:     if (xe==M*dof) {n2  = n5  = n8  = n11 = n14 = n17 = n20 = n23 = n26 = -2;}
753:   }

755:   /* If not Y periodic */
756:   if ((wrap != DA_YPERIODIC)  && (wrap != DA_XYPERIODIC) &&
757:       (wrap != DA_YZPERIODIC) && (wrap != DA_XYZPERIODIC)) {
758:     if (ys==0)   {n0  = n1  = n2  = n9  = n10 = n11 = n18 = n19 = n20 = -2;}
759:     if (ye==N)   {n6  = n7  = n8  = n15 = n16 = n17 = n24 = n25 = n26 = -2;}
760:   }

762:   /* If not Z periodic */
763:   if ((wrap != DA_ZPERIODIC)  && (wrap != DA_XZPERIODIC) &&
764:       (wrap != DA_YZPERIODIC) && (wrap != DA_XYZPERIODIC)) {
765:     if (zs==0)   {n0  = n1  = n2  = n3  = n4  = n5  = n6  = n7  = n8  = -2;}
766:     if (ze==P)   {n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;}
767:   }

769:   /* If star stencil then delete the corner neighbors */
770:   if (stencil_type == DA_STENCIL_STAR) {
771:      /* save information about corner neighbors */
772:      sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7;
773:      sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18;
774:      sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25;
775:      sn26 = n26;
776:      n0  = n1  = n2  = n3  = n5  = n6  = n7  = n8  = n9  = n11 =
777:      n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1;
778:   }


781:   PetscMalloc((Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt),&idx);
782:   PetscLogObjectMemory(da,(Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt));

784:   nn = 0;

786:   /* Bottom Level */
787:   for (k=0; k<s_z; k++) {
788:     for (i=1; i<=s_y; i++) {
789:       if (n0 >= 0) { /* left below */
790:         x_t = lx[n0 % m]*dof;
791:         y_t = ly[(n0 % (m*n))/m];
792:         z_t = lz[n0 / (m*n)];
793:         s_t = bases[n0] + x_t*y_t*z_t - (s_y-i)*x_t - s_x - (s_z-k-1)*x_t*y_t;
794:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
795:       }
796:       if (n1 >= 0) { /* directly below */
797:         x_t = x;
798:         y_t = ly[(n1 % (m*n))/m];
799:         z_t = lz[n1 / (m*n)];
800:         s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
801:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
802:       }
803:       if (n2 >= 0) { /* right below */
804:         x_t = lx[n2 % m]*dof;
805:         y_t = ly[(n2 % (m*n))/m];
806:         z_t = lz[n2 / (m*n)];
807:         s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
808:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
809:       }
810:     }

812:     for (i=0; i<y; i++) {
813:       if (n3 >= 0) { /* directly left */
814:         x_t = lx[n3 % m]*dof;
815:         y_t = y;
816:         z_t = lz[n3 / (m*n)];
817:         s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
818:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
819:       }

821:       if (n4 >= 0) { /* middle */
822:         x_t = x;
823:         y_t = y;
824:         z_t = lz[n4 / (m*n)];
825:         s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
826:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
827:       }

829:       if (n5 >= 0) { /* directly right */
830:         x_t = lx[n5 % m]*dof;
831:         y_t = y;
832:         z_t = lz[n5 / (m*n)];
833:         s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
834:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
835:       }
836:     }

838:     for (i=1; i<=s_y; i++) {
839:       if (n6 >= 0) { /* left above */
840:         x_t = lx[n6 % m]*dof;
841:         y_t = ly[(n6 % (m*n))/m];
842:         z_t = lz[n6 / (m*n)];
843:         s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
844:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
845:       }
846:       if (n7 >= 0) { /* directly above */
847:         x_t = x;
848:         y_t = ly[(n7 % (m*n))/m];
849:         z_t = lz[n7 / (m*n)];
850:         s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
851:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
852:       }
853:       if (n8 >= 0) { /* right above */
854:         x_t = lx[n8 % m]*dof;
855:         y_t = ly[(n8 % (m*n))/m];
856:         z_t = lz[n8 / (m*n)];
857:         s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
858:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
859:       }
860:     }
861:   }

863:   /* Middle Level */
864:   for (k=0; k<z; k++) {
865:     for (i=1; i<=s_y; i++) {
866:       if (n9 >= 0) { /* left below */
867:         x_t = lx[n9 % m]*dof;
868:         y_t = ly[(n9 % (m*n))/m];
869:         /* z_t = z; */
870:         s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
871:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
872:       }
873:       if (n10 >= 0) { /* directly below */
874:         x_t = x;
875:         y_t = ly[(n10 % (m*n))/m];
876:         /* z_t = z; */
877:         s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
878:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
879:       }
880:       if (n11 >= 0) { /* right below */
881:         x_t = lx[n11 % m]*dof;
882:         y_t = ly[(n11 % (m*n))/m];
883:         /* z_t = z; */
884:         s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
885:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
886:       }
887:     }

889:     for (i=0; i<y; i++) {
890:       if (n12 >= 0) { /* directly left */
891:         x_t = lx[n12 % m]*dof;
892:         y_t = y;
893:         /* z_t = z; */
894:         s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
895:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
896:       }

898:       /* Interior */
899:       s_t = bases[rank] + i*x + k*x*y;
900:       for (j=0; j<x; j++) { idx[nn++] = s_t++;}

902:       if (n14 >= 0) { /* directly right */
903:         x_t = lx[n14 % m]*dof;
904:         y_t = y;
905:         /* z_t = z; */
906:         s_t = bases[n14] + i*x_t + k*x_t*y_t;
907:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
908:       }
909:     }

911:     for (i=1; i<=s_y; i++) {
912:       if (n15 >= 0) { /* left above */
913:         x_t = lx[n15 % m]*dof;
914:         y_t = ly[(n15 % (m*n))/m];
915:         /* z_t = z; */
916:         s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
917:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
918:       }
919:       if (n16 >= 0) { /* directly above */
920:         x_t = x;
921:         y_t = ly[(n16 % (m*n))/m];
922:         /* z_t = z; */
923:         s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
924:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
925:       }
926:       if (n17 >= 0) { /* right above */
927:         x_t = lx[n17 % m]*dof;
928:         y_t = ly[(n17 % (m*n))/m];
929:         /* z_t = z; */
930:         s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
931:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
932:       }
933:     }
934:   }
935: 
936:   /* Upper Level */
937:   for (k=0; k<s_z; k++) {
938:     for (i=1; i<=s_y; i++) {
939:       if (n18 >= 0) { /* left below */
940:         x_t = lx[n18 % m]*dof;
941:         y_t = ly[(n18 % (m*n))/m];
942:         /* z_t = lz[n18 / (m*n)]; */
943:         s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
944:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
945:       }
946:       if (n19 >= 0) { /* directly below */
947:         x_t = x;
948:         y_t = ly[(n19 % (m*n))/m];
949:         /* z_t = lz[n19 / (m*n)]; */
950:         s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
951:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
952:       }
953:       if (n20 >= 0) { /* right below */
954:         x_t = lx[n20 % m]*dof;
955:         y_t = ly[(n20 % (m*n))/m];
956:         /* z_t = lz[n20 / (m*n)]; */
957:         s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
958:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
959:       }
960:     }

962:     for (i=0; i<y; i++) {
963:       if (n21 >= 0) { /* directly left */
964:         x_t = lx[n21 % m]*dof;
965:         y_t = y;
966:         /* z_t = lz[n21 / (m*n)]; */
967:         s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
968:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
969:       }

971:       if (n22 >= 0) { /* middle */
972:         x_t = x;
973:         y_t = y;
974:         /* z_t = lz[n22 / (m*n)]; */
975:         s_t = bases[n22] + i*x_t + k*x_t*y_t;
976:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
977:       }

979:       if (n23 >= 0) { /* directly right */
980:         x_t = lx[n23 % m]*dof;
981:         y_t = y;
982:         /* z_t = lz[n23 / (m*n)]; */
983:         s_t = bases[n23] + i*x_t + k*x_t*y_t;
984:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
985:       }
986:     }

988:     for (i=1; i<=s_y; i++) {
989:       if (n24 >= 0) { /* left above */
990:         x_t = lx[n24 % m]*dof;
991:         y_t = ly[(n24 % (m*n))/m];
992:         /* z_t = lz[n24 / (m*n)]; */
993:         s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
994:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
995:       }
996:       if (n25 >= 0) { /* directly above */
997:         x_t = x;
998:         y_t = ly[(n25 % (m*n))/m];
999:         /* z_t = lz[n25 / (m*n)]; */
1000:         s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
1001:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1002:       }
1003:       if (n26 >= 0) { /* right above */
1004:         x_t = lx[n26 % m]*dof;
1005:         y_t = ly[(n26 % (m*n))/m];
1006:         /* z_t = lz[n26 / (m*n)]; */
1007:         s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
1008:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1009:       }
1010:     }
1011:   }
1012:   base = bases[rank];
1013:   ISCreateGeneral(comm,nn,idx,&from);
1014:   VecScatterCreate(global,from,local,to,&gtol);
1015:   PetscLogObjectParent(da,gtol);
1016:   PetscLogObjectParent(da,to);
1017:   PetscLogObjectParent(da,from);
1018:   ISDestroy(to);
1019:   ISDestroy(from);
1020:   da->stencil_type = stencil_type;
1021:   da->M  = M;  da->N  = N; da->P = P;
1022:   da->m  = m;  da->n  = n; da->p = p;
1023:   da->w  = dof;  da->s  = s;
1024:   da->xs = xs; da->xe = xe; da->ys = ys; da->ye = ye; da->zs = zs; da->ze = ze;
1025:   da->Xs = Xs; da->Xe = Xe; da->Ys = Ys; da->Ye = Ye; da->Zs = Zs; da->Ze = Ze;

1027:   VecDestroy(local);
1028:   VecDestroy(global);

1030:   if (stencil_type == DA_STENCIL_STAR) {
1031:     /*
1032:         Recompute the local to global mappings, this time keeping the 
1033:       information about the cross corner processor numbers.
1034:     */
1035:     n0  = sn0;  n1  = sn1;  n2  = sn2;  n3  = sn3;  n5  = sn5;  n6  = sn6; n7 = sn7;
1036:     n8  = sn8;  n9  = sn9;  n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18;
1037:     n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25;
1038:     n26 = sn26;

1040:     nn = 0;

1042:     /* Bottom Level */
1043:     for (k=0; k<s_z; k++) {
1044:       for (i=1; i<=s_y; i++) {
1045:         if (n0 >= 0) { /* left below */
1046:           x_t = lx[n0 % m]*dof;
1047:           y_t = ly[(n0 % (m*n))/m];
1048:           z_t = lz[n0 / (m*n)];
1049:           s_t = bases[n0] + x_t*y_t*z_t - (s_y-i)*x_t - s_x - (s_z-k-1)*x_t*y_t;
1050:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1051:         }
1052:         if (n1 >= 0) { /* directly below */
1053:           x_t = x;
1054:           y_t = ly[(n1 % (m*n))/m];
1055:           z_t = lz[n1 / (m*n)];
1056:           s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
1057:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1058:         }
1059:         if (n2 >= 0) { /* right below */
1060:           x_t = lx[n2 % m]*dof;
1061:           y_t = ly[(n2 % (m*n))/m];
1062:           z_t = lz[n2 / (m*n)];
1063:           s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
1064:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1065:         }
1066:       }

1068:       for (i=0; i<y; i++) {
1069:         if (n3 >= 0) { /* directly left */
1070:           x_t = lx[n3 % m]*dof;
1071:           y_t = y;
1072:           z_t = lz[n3 / (m*n)];
1073:           s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1074:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1075:         }

1077:         if (n4 >= 0) { /* middle */
1078:           x_t = x;
1079:           y_t = y;
1080:           z_t = lz[n4 / (m*n)];
1081:           s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1082:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1083:         }

1085:         if (n5 >= 0) { /* directly right */
1086:           x_t = lx[n5 % m]*dof;
1087:           y_t = y;
1088:           z_t = lz[n5 / (m*n)];
1089:           s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1090:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1091:         }
1092:       }

1094:       for (i=1; i<=s_y; i++) {
1095:         if (n6 >= 0) { /* left above */
1096:           x_t = lx[n6 % m]*dof;
1097:           y_t = ly[(n6 % (m*n))/m];
1098:           z_t = lz[n6 / (m*n)];
1099:           s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1100:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1101:         }
1102:         if (n7 >= 0) { /* directly above */
1103:           x_t = x;
1104:           y_t = ly[(n7 % (m*n))/m];
1105:           z_t = lz[n7 / (m*n)];
1106:           s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1107:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1108:         }
1109:         if (n8 >= 0) { /* right above */
1110:           x_t = lx[n8 % m]*dof;
1111:           y_t = ly[(n8 % (m*n))/m];
1112:           z_t = lz[n8 / (m*n)];
1113:           s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1114:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1115:         }
1116:       }
1117:     }

1119:     /* Middle Level */
1120:     for (k=0; k<z; k++) {
1121:       for (i=1; i<=s_y; i++) {
1122:         if (n9 >= 0) { /* left below */
1123:           x_t = lx[n9 % m]*dof;
1124:           y_t = ly[(n9 % (m*n))/m];
1125:           /* z_t = z; */
1126:           s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
1127:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1128:         }
1129:         if (n10 >= 0) { /* directly below */
1130:           x_t = x;
1131:           y_t = ly[(n10 % (m*n))/m];
1132:           /* z_t = z; */
1133:           s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1134:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1135:         }
1136:         if (n11 >= 0) { /* right below */
1137:           x_t = lx[n11 % m]*dof;
1138:           y_t = ly[(n11 % (m*n))/m];
1139:           /* z_t = z; */
1140:           s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1141:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1142:         }
1143:       }

1145:       for (i=0; i<y; i++) {
1146:         if (n12 >= 0) { /* directly left */
1147:           x_t = lx[n12 % m]*dof;
1148:           y_t = y;
1149:           /* z_t = z; */
1150:           s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
1151:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1152:         }

1154:         /* Interior */
1155:         s_t = bases[rank] + i*x + k*x*y;
1156:         for (j=0; j<x; j++) { idx[nn++] = s_t++;}

1158:         if (n14 >= 0) { /* directly right */
1159:           x_t = lx[n14 % m]*dof;
1160:           y_t = y;
1161:           /* z_t = z; */
1162:           s_t = bases[n14] + i*x_t + k*x_t*y_t;
1163:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1164:         }
1165:       }

1167:       for (i=1; i<=s_y; i++) {
1168:         if (n15 >= 0) { /* left above */
1169:           x_t = lx[n15 % m]*dof;
1170:           y_t = ly[(n15 % (m*n))/m];
1171:           /* z_t = z; */
1172:           s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
1173:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1174:         }
1175:         if (n16 >= 0) { /* directly above */
1176:           x_t = x;
1177:           y_t = ly[(n16 % (m*n))/m];
1178:           /* z_t = z; */
1179:           s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
1180:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1181:         }
1182:         if (n17 >= 0) { /* right above */
1183:           x_t = lx[n17 % m]*dof;
1184:           y_t = ly[(n17 % (m*n))/m];
1185:           /* z_t = z; */
1186:           s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
1187:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1188:         }
1189:       }
1190:     }
1191: 
1192:     /* Upper Level */
1193:     for (k=0; k<s_z; k++) {
1194:       for (i=1; i<=s_y; i++) {
1195:         if (n18 >= 0) { /* left below */
1196:           x_t = lx[n18 % m]*dof;
1197:           y_t = ly[(n18 % (m*n))/m];
1198:           /* z_t = lz[n18 / (m*n)]; */
1199:           s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
1200:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1201:         }
1202:         if (n19 >= 0) { /* directly below */
1203:           x_t = x;
1204:           y_t = ly[(n19 % (m*n))/m];
1205:           /* z_t = lz[n19 / (m*n)]; */
1206:           s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1207:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1208:         }
1209:         if (n20 >= 0) { /* right below */
1210:           x_t = lx[n20 % m]*dof;
1211:           y_t = ly[(n20 % (m*n))/m];
1212:           /* z_t = lz[n20 / (m*n)]; */
1213:           s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1214:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1215:         }
1216:       }

1218:       for (i=0; i<y; i++) {
1219:         if (n21 >= 0) { /* directly left */
1220:           x_t = lx[n21 % m]*dof;
1221:           y_t = y;
1222:           /* z_t = lz[n21 / (m*n)]; */
1223:           s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
1224:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1225:         }

1227:         if (n22 >= 0) { /* middle */
1228:           x_t = x;
1229:           y_t = y;
1230:           /* z_t = lz[n22 / (m*n)]; */
1231:           s_t = bases[n22] + i*x_t + k*x_t*y_t;
1232:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1233:         }

1235:         if (n23 >= 0) { /* directly right */
1236:           x_t = lx[n23 % m]*dof;
1237:           y_t = y;
1238:           /* z_t = lz[n23 / (m*n)]; */
1239:           s_t = bases[n23] + i*x_t + k*x_t*y_t;
1240:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1241:         }
1242:       }

1244:       for (i=1; i<=s_y; i++) {
1245:         if (n24 >= 0) { /* left above */
1246:           x_t = lx[n24 % m]*dof;
1247:           y_t = ly[(n24 % (m*n))/m];
1248:           /* z_t = lz[n24 / (m*n)]; */
1249:           s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
1250:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1251:         }
1252:         if (n25 >= 0) { /* directly above */
1253:           x_t = x;
1254:           y_t = ly[(n25 % (m*n))/m];
1255:           /* z_t = lz[n25 / (m*n)]; */
1256:           s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
1257:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1258:         }
1259:         if (n26 >= 0) { /* right above */
1260:           x_t = lx[n26 % m]*dof;
1261:           y_t = ly[(n26 % (m*n))/m];
1262:           /* z_t = lz[n26 / (m*n)]; */
1263:           s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
1264:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1265:         }
1266:       }
1267:     }
1268:   }
1269:   /* redo idx to include "missing" ghost points */
1270:   /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
1271: 
1272:   /* Assume Nodes are Internal to the Cube */
1273: 
1274:   n0  = rank - m*n - m - 1;
1275:   n1  = rank - m*n - m;
1276:   n2  = rank - m*n - m + 1;
1277:   n3  = rank - m*n -1;
1278:   n4  = rank - m*n;
1279:   n5  = rank - m*n + 1;
1280:   n6  = rank - m*n + m - 1;
1281:   n7  = rank - m*n + m;
1282:   n8  = rank - m*n + m + 1;

1284:   n9  = rank - m - 1;
1285:   n10 = rank - m;
1286:   n11 = rank - m + 1;
1287:   n12 = rank - 1;
1288:   n14 = rank + 1;
1289:   n15 = rank + m - 1;
1290:   n16 = rank + m;
1291:   n17 = rank + m + 1;

1293:   n18 = rank + m*n - m - 1;
1294:   n19 = rank + m*n - m;
1295:   n20 = rank + m*n - m + 1;
1296:   n21 = rank + m*n - 1;
1297:   n22 = rank + m*n;
1298:   n23 = rank + m*n + 1;
1299:   n24 = rank + m*n + m - 1;
1300:   n25 = rank + m*n + m;
1301:   n26 = rank + m*n + m + 1;

1303:   /* Assume Pieces are on Faces of Cube */

1305:   if (xs == 0) { /* First assume not corner or edge */
1306:     n0  = rank       -1 - (m*n);
1307:     n3  = rank + m   -1 - (m*n);
1308:     n6  = rank + 2*m -1 - (m*n);
1309:     n9  = rank       -1;
1310:     n12 = rank + m   -1;
1311:     n15 = rank + 2*m -1;
1312:     n18 = rank       -1 + (m*n);
1313:     n21 = rank + m   -1 + (m*n);
1314:     n24 = rank + 2*m -1 + (m*n);
1315:    }

1317:   if (xe == M*dof) { /* First assume not corner or edge */
1318:     n2  = rank -2*m +1 - (m*n);
1319:     n5  = rank - m  +1 - (m*n);
1320:     n8  = rank      +1 - (m*n);
1321:     n11 = rank -2*m +1;
1322:     n14 = rank - m  +1;
1323:     n17 = rank      +1;
1324:     n20 = rank -2*m +1 + (m*n);
1325:     n23 = rank - m  +1 + (m*n);
1326:     n26 = rank      +1 + (m*n);
1327:   }

1329:   if (ys==0) { /* First assume not corner or edge */
1330:     n0  = rank + m * (n-1) -1 - (m*n);
1331:     n1  = rank + m * (n-1)    - (m*n);
1332:     n2  = rank + m * (n-1) +1 - (m*n);
1333:     n9  = rank + m * (n-1) -1;
1334:     n10 = rank + m * (n-1);
1335:     n11 = rank + m * (n-1) +1;
1336:     n18 = rank + m * (n-1) -1 + (m*n);
1337:     n19 = rank + m * (n-1)    + (m*n);
1338:     n20 = rank + m * (n-1) +1 + (m*n);
1339:   }

1341:   if (ye == N) { /* First assume not corner or edge */
1342:     n6  = rank - m * (n-1) -1 - (m*n);
1343:     n7  = rank - m * (n-1)    - (m*n);
1344:     n8  = rank - m * (n-1) +1 - (m*n);
1345:     n15 = rank - m * (n-1) -1;
1346:     n16 = rank - m * (n-1);
1347:     n17 = rank - m * (n-1) +1;
1348:     n24 = rank - m * (n-1) -1 + (m*n);
1349:     n25 = rank - m * (n-1)    + (m*n);
1350:     n26 = rank - m * (n-1) +1 + (m*n);
1351:   }
1352: 
1353:   if (zs == 0) { /* First assume not corner or edge */
1354:     n0 = size - (m*n) + rank - m - 1;
1355:     n1 = size - (m*n) + rank - m;
1356:     n2 = size - (m*n) + rank - m + 1;
1357:     n3 = size - (m*n) + rank - 1;
1358:     n4 = size - (m*n) + rank;
1359:     n5 = size - (m*n) + rank + 1;
1360:     n6 = size - (m*n) + rank + m - 1;
1361:     n7 = size - (m*n) + rank + m ;
1362:     n8 = size - (m*n) + rank + m + 1;
1363:   }

1365:   if (ze == P) { /* First assume not corner or edge */
1366:     n18 = (m*n) - (size-rank) - m - 1;
1367:     n19 = (m*n) - (size-rank) - m;
1368:     n20 = (m*n) - (size-rank) - m + 1;
1369:     n21 = (m*n) - (size-rank) - 1;
1370:     n22 = (m*n) - (size-rank);
1371:     n23 = (m*n) - (size-rank) + 1;
1372:     n24 = (m*n) - (size-rank) + m - 1;
1373:     n25 = (m*n) - (size-rank) + m;
1374:     n26 = (m*n) - (size-rank) + m + 1;
1375:   }

1377:   if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */
1378:     n0 = size - m*n + rank + m-1 - m;
1379:     n3 = size - m*n + rank + m-1;
1380:     n6 = size - m*n + rank + m-1 + m;
1381:   }
1382: 
1383:   if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */
1384:     n18 = m*n - (size - rank) + m-1 - m;
1385:     n21 = m*n - (size - rank) + m-1;
1386:     n24 = m*n - (size - rank) + m-1 + m;
1387:   }

1389:   if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */
1390:     n0  = rank + m*n -1 - m*n;
1391:     n9  = rank + m*n -1;
1392:     n18 = rank + m*n -1 + m*n;
1393:   }

1395:   if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */
1396:     n6  = rank - m*(n-1) + m-1 - m*n;
1397:     n15 = rank - m*(n-1) + m-1;
1398:     n24 = rank - m*(n-1) + m-1 + m*n;
1399:   }

1401:   if ((xe==M*dof) && (zs==0)) { /* Assume an edge, not corner */
1402:     n2 = size - (m*n-rank) - (m-1) - m;
1403:     n5 = size - (m*n-rank) - (m-1);
1404:     n8 = size - (m*n-rank) - (m-1) + m;
1405:   }

1407:   if ((xe==M*dof) && (ze==P)) { /* Assume an edge, not corner */
1408:     n20 = m*n - (size - rank) - (m-1) - m;
1409:     n23 = m*n - (size - rank) - (m-1);
1410:     n26 = m*n - (size - rank) - (m-1) + m;
1411:   }

1413:   if ((xe==M*dof) && (ys==0)) { /* Assume an edge, not corner */
1414:     n2  = rank + m*(n-1) - (m-1) - m*n;
1415:     n11 = rank + m*(n-1) - (m-1);
1416:     n20 = rank + m*(n-1) - (m-1) + m*n;
1417:   }

1419:   if ((xe==M*dof) && (ye==N)) { /* Assume an edge, not corner */
1420:     n8  = rank - m*n +1 - m*n;
1421:     n17 = rank - m*n +1;
1422:     n26 = rank - m*n +1 + m*n;
1423:   }

1425:   if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */
1426:     n0 = size - m + rank -1;
1427:     n1 = size - m + rank;
1428:     n2 = size - m + rank +1;
1429:   }

1431:   if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */
1432:     n18 = m*n - (size - rank) + m*(n-1) -1;
1433:     n19 = m*n - (size - rank) + m*(n-1);
1434:     n20 = m*n - (size - rank) + m*(n-1) +1;
1435:   }

1437:   if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */
1438:     n6 = size - (m*n-rank) - m * (n-1) -1;
1439:     n7 = size - (m*n-rank) - m * (n-1);
1440:     n8 = size - (m*n-rank) - m * (n-1) +1;
1441:   }

1443:   if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */
1444:     n24 = rank - (size-m) -1;
1445:     n25 = rank - (size-m);
1446:     n26 = rank - (size-m) +1;
1447:   }

1449:   /* Check for Corners */
1450:   if ((xs==0)   && (ys==0) && (zs==0)) { n0  = size -1;}
1451:   if ((xs==0)   && (ys==0) && (ze==P)) { n18 = m*n-1;}
1452:   if ((xs==0)   && (ye==N) && (zs==0)) { n6  = (size-1)-m*(n-1);}
1453:   if ((xs==0)   && (ye==N) && (ze==P)) { n24 = m-1;}
1454:   if ((xe==M*dof) && (ys==0) && (zs==0)) { n2  = size-m;}
1455:   if ((xe==M*dof) && (ys==0) && (ze==P)) { n20 = m*n-m;}
1456:   if ((xe==M*dof) && (ye==N) && (zs==0)) { n8  = size-m*n;}
1457:   if ((xe==M*dof) && (ye==N) && (ze==P)) { n26 = 0;}

1459:   /* Check for when not X,Y, and Z Periodic */

1461:   /* If not X periodic */
1462:   if (!DAXPeriodic(wrap)){
1463:     if (xs==0)   {n0  = n3  = n6  = n9  = n12 = n15 = n18 = n21 = n24 = -2;}
1464:     if (xe==M*dof) {n2  = n5  = n8  = n11 = n14 = n17 = n20 = n23 = n26 = -2;}
1465:   }

1467:   /* If not Y periodic */
1468:   if (!DAYPeriodic(wrap)){
1469:     if (ys==0)   {n0  = n1  = n2  = n9  = n10 = n11 = n18 = n19 = n20 = -2;}
1470:     if (ye==N)   {n6  = n7  = n8  = n15 = n16 = n17 = n24 = n25 = n26 = -2;}
1471:   }

1473:   /* If not Z periodic */
1474:   if (!DAZPeriodic(wrap)){
1475:     if (zs==0)   {n0  = n1  = n2  = n3  = n4  = n5  = n6  = n7  = n8  = -2;}
1476:     if (ze==P)   {n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;}
1477:   }

1479:   nn = 0;

1481:   /* Bottom Level */
1482:   for (k=0; k<s_z; k++) {
1483:     for (i=1; i<=s_y; i++) {
1484:       if (n0 >= 0) { /* left below */
1485:         x_t = lx[n0 % m]*dof;
1486:         y_t = ly[(n0 % (m*n))/m];
1487:         z_t = lz[n0 / (m*n)];
1488:         s_t = bases[n0] + x_t*y_t*z_t - (s_y-i)*x_t -s_x - (s_z-k-1)*x_t*y_t;
1489:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1490:       }
1491:       if (n1 >= 0) { /* directly below */
1492:         x_t = x;
1493:         y_t = ly[(n1 % (m*n))/m];
1494:         z_t = lz[n1 / (m*n)];
1495:         s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
1496:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1497:       }
1498:       if (n2 >= 0) { /* right below */
1499:         x_t = lx[n2 % m]*dof;
1500:         y_t = ly[(n2 % (m*n))/m];
1501:         z_t = lz[n2 / (m*n)];
1502:         s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
1503:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1504:       }
1505:     }

1507:     for (i=0; i<y; i++) {
1508:       if (n3 >= 0) { /* directly left */
1509:         x_t = lx[n3 % m]*dof;
1510:         y_t = y;
1511:         z_t = lz[n3 / (m*n)];
1512:         s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1513:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1514:       }

1516:       if (n4 >= 0) { /* middle */
1517:         x_t = x;
1518:         y_t = y;
1519:         z_t = lz[n4 / (m*n)];
1520:         s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1521:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1522:       }

1524:       if (n5 >= 0) { /* directly right */
1525:         x_t = lx[n5 % m]*dof;
1526:         y_t = y;
1527:         z_t = lz[n5 / (m*n)];
1528:         s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1529:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1530:       }
1531:     }

1533:     for (i=1; i<=s_y; i++) {
1534:       if (n6 >= 0) { /* left above */
1535:         x_t = lx[n6 % m]*dof;
1536:         y_t = ly[(n6 % (m*n))/m];
1537:         z_t = lz[n6 / (m*n)];
1538:         s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1539:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1540:       }
1541:       if (n7 >= 0) { /* directly above */
1542:         x_t = x;
1543:         y_t = ly[(n7 % (m*n))/m];
1544:         z_t = lz[n7 / (m*n)];
1545:         s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1546:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1547:       }
1548:       if (n8 >= 0) { /* right above */
1549:         x_t = lx[n8 % m]*dof;
1550:         y_t = ly[(n8 % (m*n))/m];
1551:         z_t = lz[n8 / (m*n)];
1552:         s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1553:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1554:       }
1555:     }
1556:   }

1558:   /* Middle Level */
1559:   for (k=0; k<z; k++) {
1560:     for (i=1; i<=s_y; i++) {
1561:       if (n9 >= 0) { /* left below */
1562:         x_t = lx[n9 % m]*dof;
1563:         y_t = ly[(n9 % (m*n))/m];
1564:         /* z_t = z; */
1565:         s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
1566:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1567:       }
1568:       if (n10 >= 0) { /* directly below */
1569:         x_t = x;
1570:         y_t = ly[(n10 % (m*n))/m];
1571:         /* z_t = z; */
1572:         s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1573:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1574:       }
1575:       if (n11 >= 0) { /* right below */
1576:         x_t = lx[n11 % m]*dof;
1577:         y_t = ly[(n11 % (m*n))/m];
1578:         /* z_t = z; */
1579:         s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1580:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1581:       }
1582:     }

1584:     for (i=0; i<y; i++) {
1585:       if (n12 >= 0) { /* directly left */
1586:         x_t = lx[n12 % m]*dof;
1587:         y_t = y;
1588:         /* z_t = z; */
1589:         s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
1590:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1591:       }

1593:       /* Interior */
1594:       s_t = bases[rank] + i*x + k*x*y;
1595:       for (j=0; j<x; j++) { idx[nn++] = s_t++;}

1597:       if (n14 >= 0) { /* directly right */
1598:         x_t = lx[n14 % m]*dof;
1599:         y_t = y;
1600:         /* z_t = z; */
1601:         s_t = bases[n14] + i*x_t + k*x_t*y_t;
1602:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1603:       }
1604:     }

1606:     for (i=1; i<=s_y; i++) {
1607:       if (n15 >= 0) { /* left above */
1608:         x_t = lx[n15 % m]*dof;
1609:         y_t = ly[(n15 % (m*n))/m];
1610:         /* z_t = z; */
1611:         s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
1612:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1613:       }
1614:       if (n16 >= 0) { /* directly above */
1615:         x_t = x;
1616:         y_t = ly[(n16 % (m*n))/m];
1617:         /* z_t = z; */
1618:         s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
1619:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1620:       }
1621:       if (n17 >= 0) { /* right above */
1622:         x_t = lx[n17 % m]*dof;
1623:         y_t = ly[(n17 % (m*n))/m];
1624:         /* z_t = z; */
1625:         s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
1626:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1627:       }
1628:     }
1629:   }
1630: 
1631:   /* Upper Level */
1632:   for (k=0; k<s_z; k++) {
1633:     for (i=1; i<=s_y; i++) {
1634:       if (n18 >= 0) { /* left below */
1635:         x_t = lx[n18 % m]*dof;
1636:         y_t = ly[(n18 % (m*n))/m];
1637:         /* z_t = lz[n18 / (m*n)]; */
1638:         s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
1639:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1640:       }
1641:       if (n19 >= 0) { /* directly below */
1642:         x_t = x;
1643:         y_t = ly[(n19 % (m*n))/m];
1644:         /* z_t = lz[n19 / (m*n)]; */
1645:         s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1646:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1647:       }
1648:       if (n20 >= 0) { /* right belodof */
1649:         x_t = lx[n20 % m]*dof;
1650:         y_t = ly[(n20 % (m*n))/m];
1651:         /* z_t = lz[n20 / (m*n)]; */
1652:         s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1653:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1654:       }
1655:     }

1657:     for (i=0; i<y; i++) {
1658:       if (n21 >= 0) { /* directly left */
1659:         x_t = lx[n21 % m]*dof;
1660:         y_t = y;
1661:         /* z_t = lz[n21 / (m*n)]; */
1662:         s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
1663:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1664:       }

1666:       if (n22 >= 0) { /* middle */
1667:         x_t = x;
1668:         y_t = y;
1669:         /* z_t = lz[n22 / (m*n)]; */
1670:         s_t = bases[n22] + i*x_t + k*x_t*y_t;
1671:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1672:       }

1674:       if (n23 >= 0) { /* directly right */
1675:         x_t = lx[n23 % m]*dof;
1676:         y_t = y;
1677:         /* z_t = lz[n23 / (m*n)]; */
1678:         s_t = bases[n23] + i*x_t + k*x_t*y_t;
1679:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1680:       }
1681:     }

1683:     for (i=1; i<=s_y; i++) {
1684:       if (n24 >= 0) { /* left above */
1685:         x_t = lx[n24 % m]*dof;
1686:         y_t = ly[(n24 % (m*n))/m];
1687:         /* z_t = lz[n24 / (m*n)]; */
1688:         s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
1689:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1690:       }
1691:       if (n25 >= 0) { /* directly above */
1692:         x_t = x;
1693:         y_t = ly[(n25 % (m*n))/m];
1694:         /* z_t = lz[n25 / (m*n)]; */
1695:         s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
1696:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1697:       }
1698:       if (n26 >= 0) { /* right above */
1699:         x_t = lx[n26 % m]*dof;
1700:         y_t = ly[(n26 % (m*n))/m];
1701:         /* z_t = lz[n26 / (m*n)]; */
1702:         s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
1703:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1704:       }
1705:     }
1706:   }
1707:   PetscFree(bases);
1708:   da->gtol      = gtol;
1709:   da->ltog      = ltog;
1710:   da->idx       = idx;
1711:   da->Nl        = nn;
1712:   da->base      = base;
1713:   da->ops->view = DAView_3d;
1714:   da->wrap      = wrap;
1715:   *inra = da;

1717:   /* 
1718:      Set the local to global ordering in the global vector, this allows use
1719:      of VecSetValuesLocal().
1720:   */
1721:   ISLocalToGlobalMappingCreateNC(comm,nn,idx,&da->ltogmap);
1722:   ISLocalToGlobalMappingBlock(da->ltogmap,da->w,&da->ltogmapb);
1723:   PetscLogObjectParent(da,da->ltogmap);

1725:   da->ltol = PETSC_NULL;
1726:   da->ao   = PETSC_NULL;

1728:   if (!flx) {
1729:     PetscMalloc(m*sizeof(PetscInt),&flx);
1730:     PetscMemcpy(flx,lx,m*sizeof(PetscInt));
1731:   }
1732:   if (!fly) {
1733:     PetscMalloc(n*sizeof(PetscInt),&fly);
1734:     PetscMemcpy(fly,ly,n*sizeof(PetscInt));
1735:   }
1736:   if (!flz) {
1737:     PetscMalloc(p*sizeof(PetscInt),&flz);
1738:     PetscMemcpy(flz,lz,p*sizeof(PetscInt));
1739:   }
1740:   da->lx = flx;
1741:   da->ly = fly;
1742:   da->lz = flz;

1744:   PetscOptionsHasName(PETSC_NULL,"-da_view",&flg1);
1745:   if (flg1) {DAView(da,PETSC_VIEWER_STDOUT_(da->comm));}
1746:   PetscOptionsHasName(PETSC_NULL,"-da_view_draw",&flg1);
1747:   if (flg1) {DAView(da,PETSC_VIEWER_DRAW_(da->comm));}
1748:   PetscOptionsHasName(PETSC_NULL,"-help",&flg1);
1749:   if (flg1) {DAPrintHelp(da);}
1750:   PetscPublishAll(da);

1752:   return(0);
1753: }