Actual source code: da3.c

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

 7:  #include private/daimpl.h

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

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

 20:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
 21:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
 22:   if (iascii) {
 23:     PetscViewerFormat format;

 25:     PetscViewerGetFormat(viewer, &format);
 26:     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) {
 27:       DALocalInfo info;
 28:       DAGetLocalInfo(da,&info);
 29:       PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %D N %D P %D m %D n %D p %D w %D s %D\n",
 30:                                                 rank,da->M,da->N,da->P,da->m,da->n,da->p,da->w,da->s);
 31:       PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D, Y range of indices: %D %D, Z range of indices: %D %D\n",
 32:                                                 info.xs,info.xs+info.xm,info.ys,info.ys+info.ym,info.zs,info.zs+info.zm);
 33: #if !defined(PETSC_USE_COMPLEX)
 34:       if (da->coordinates) {
 35:         PetscInt  last;
 36:         PetscReal *coors;
 37:         VecGetArray(da->coordinates,&coors);
 38:         VecGetLocalSize(da->coordinates,&last);
 39:         last = last - 3;
 40:         PetscViewerASCIISynchronizedPrintf(viewer,"Lower left corner %G %G %G : Upper right %G %G %G\n",
 41:                                                   coors[0],coors[1],coors[2],coors[last],coors[last+1],coors[last+2]);
 42:         VecRestoreArray(da->coordinates,&coors);
 43:       }
 44: #endif
 45:       PetscViewerFlush(viewer);
 46:     }
 47:   } else if (isdraw) {
 48:     PetscDraw       draw;
 49:     PetscReal     ymin = -1.0,ymax = (PetscReal)da->N;
 50:     PetscReal     xmin = -1.0,xmax = (PetscReal)((da->M+2)*da->P),x,y,ycoord,xcoord;
 51:     PetscInt        k,plane,base,*idx;
 52:     char       node[10];
 53:     PetscTruth isnull;

 55:     PetscViewerDrawGetDraw(viewer,0,&draw);
 56:     PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
 57:     PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);
 58:     PetscDrawSynchronizedClear(draw);

 60:     /* first processor draw all node lines */
 61:     if (!rank) {
 62:       for (k=0; k<da->P; k++) {
 63:         ymin = 0.0; ymax = (PetscReal)(da->N - 1);
 64:         for (xmin=(PetscReal)(k*(da->M+1)); xmin<(PetscReal)(da->M+(k*(da->M+1))); xmin++) {
 65:           PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);
 66:         }
 67: 
 68:         xmin = (PetscReal)(k*(da->M+1)); xmax = xmin + (PetscReal)(da->M - 1);
 69:         for (ymin=0; ymin<(PetscReal)da->N; ymin++) {
 70:           PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);
 71:         }
 72:       }
 73:     }
 74:     PetscDrawSynchronizedFlush(draw);
 75:     PetscDrawPause(draw);

 77:     for (k=0; k<da->P; k++) {  /*Go through and draw for each plane*/
 78:       if ((k >= da->zs) && (k < da->ze)) {
 79:         /* draw my box */
 80:         ymin = da->ys;
 81:         ymax = da->ye - 1;
 82:         xmin = da->xs/da->w    + (da->M+1)*k;
 83:         xmax =(da->xe-1)/da->w + (da->M+1)*k;

 85:         PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);
 86:         PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);
 87:         PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);
 88:         PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);

 90:         xmin = da->xs/da->w;
 91:         xmax =(da->xe-1)/da->w;

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

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

100:         for (y=ymin; y<=ymax; y++) {
101:           for (x=xmin+(da->M+1)*k; x<=xmax+(da->M+1)*k; x++) {
102:             sprintf(node,"%d",(int)base++);
103:             PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);
104:           }
105:         }
106: 
107:       }
108:     }
109:     PetscDrawSynchronizedFlush(draw);
110:     PetscDrawPause(draw);

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

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

135:             if (x<xmin)  { xcoord = xmax - (xmin-x); }
136:             if (x>=xmax) { xcoord = xmin + (x-xmax); }
137:             PetscDrawString(draw,xcoord/da->w,ycoord,PETSC_DRAW_BLUE,node);
138:             base+=da->w;
139:           }
140:         }
141:       }
142:     }
143:     PetscDrawSynchronizedFlush(draw);
144:     PetscDrawPause(draw);
145:   } else {
146:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for DA 3d",((PetscObject)viewer)->type_name);
147:   }
148:   return(0);
149: }

154: PetscErrorCode  DACreate_3D(DA da)
155: {
156:   const PetscInt       dim          = da->dim;
157:   const PetscInt       M            = da->M;
158:   const PetscInt       N            = da->N;
159:   const PetscInt       P            = da->P;
160:   PetscInt             m            = da->m;
161:   PetscInt             n            = da->n;
162:   PetscInt             p            = da->p;
163:   const PetscInt       dof          = da->w;
164:   const PetscInt       s            = da->s;
165:   const DAPeriodicType wrap         = da->wrap;
166:   const DAStencilType  stencil_type = da->stencil_type;
167:   PetscInt            *lx           = da->lx;
168:   PetscInt            *ly           = da->ly;
169:   PetscInt            *lz           = da->lz;
170:   MPI_Comm             comm;
171:   PetscMPIInt    rank,size;
172:   PetscInt       xs = 0,xe,ys = 0,ye,zs = 0,ze,x = 0,y = 0,z = 0,Xs,Xe,Ys,Ye,Zs,Ze,start,end,pm;
173:   PetscInt       left,up,down,bottom,top,i,j,k,*idx,nn;
174:   PetscInt       n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14;
175:   PetscInt       n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26;
176:   PetscInt       *bases,*ldims,x_t,y_t,z_t,s_t,base,count,s_x,s_y,s_z;
177:   PetscInt       sn0 = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0;
178:   PetscInt       sn8 = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0;
179:   PetscInt       sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0;
180:   Vec            local,global;
181:   VecScatter     ltog,gtol;
182:   IS             to,from;

186: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
187:   DMInitializePackage(PETSC_NULL);
188: #endif

190:   if (dim != PETSC_DECIDE && dim != 3) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Dimension should be 3: %D",dim);
191:   if (dof < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
192:   if (s < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);

194:   PetscObjectGetComm((PetscObject) da, &comm);
195:   MPI_Comm_size(comm,&size);
196:   MPI_Comm_rank(comm,&rank);

198:   da->dim = 3;
199:   PetscMalloc(dof*sizeof(char*),&da->fieldname);
200:   PetscMemzero(da->fieldname,dof*sizeof(char*));

202:   if (m != PETSC_DECIDE) {
203:     if (m < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);}
204:     else if (m > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);}
205:   }
206:   if (n != PETSC_DECIDE) {
207:     if (n < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);}
208:     else if (n > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);}
209:   }
210:   if (p != PETSC_DECIDE) {
211:     if (p < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %D",p);}
212:     else if (p > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %D %d",p,size);}
213:   }

215:   /* Partition the array among the processors */
216:   if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) {
217:     m = size/(n*p);
218:   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
219:     n = size/(m*p);
220:   } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
221:     p = size/(m*n);
222:   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
223:     /* try for squarish distribution */
224:     m = (int)(0.5 + sqrt(((PetscReal)M)*((PetscReal)size)/((PetscReal)N*p)));
225:     if (!m) m = 1;
226:     while (m > 0) {
227:       n = size/(m*p);
228:       if (m*n*p == size) break;
229:       m--;
230:     }
231:     if (!m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %D",p);
232:     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
233:   } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
234:     /* try for squarish distribution */
235:     m = (int)(0.5 + sqrt(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n)));
236:     if (!m) m = 1;
237:     while (m > 0) {
238:       p = size/(m*n);
239:       if (m*n*p == size) break;
240:       m--;
241:     }
242:     if (!m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %D",n);
243:     if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;}
244:   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
245:     /* try for squarish distribution */
246:     n = (int)(0.5 + sqrt(((PetscReal)N)*((PetscReal)size)/((PetscReal)P*m)));
247:     if (!n) n = 1;
248:     while (n > 0) {
249:       p = size/(m*n);
250:       if (m*n*p == size) break;
251:       n--;
252:     }
253:     if (!n) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %D",n);
254:     if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;}
255:   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
256:     /* try for squarish distribution */
257:     n = (PetscInt)(0.5 + pow(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),(PetscReal)(1./3.)));
258:     if (!n) n = 1;
259:     while (n > 0) {
260:       pm = size/n;
261:       if (n*pm == size) break;
262:       n--;
263:     }
264:     if (!n) n = 1;
265:     m = (PetscInt)(0.5 + sqrt(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n)));
266:     if (!m) m = 1;
267:     while (m > 0) {
268:       p = size/(m*n);
269:       if (m*n*p == size) break;
270:       m--;
271:     }
272:     if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;}
273:   } else if (m*n*p != size) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");

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

280:   /* 
281:      Determine locally owned region 
282:      [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes 
283:   */

285:   if (!lx) {
286:     PetscMalloc(m*sizeof(PetscInt), &da->lx);
287:     lx = da->lx;
288:     for (i=0; i<m; i++) {
289:       lx[i] = M/m + ((M % m) > (i % m));
290:     }
291:   }
292:   x  = lx[rank % m];
293:   xs = 0;
294:   for (i=0; i<(rank%m); i++) { xs += lx[i];}
295:   if (m > 1 && x < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column width is too thin for stencil! %D %D",x,s);

297:   if (!ly) {
298:     PetscMalloc(n*sizeof(PetscInt), &da->ly);
299:     ly = da->ly;
300:     for (i=0; i<n; i++) {
301:       ly[i] = N/n + ((N % n) > (i % n));
302:     }
303:   }
304:   y  = ly[(rank % (m*n))/m];
305:   if (n > 1 && y < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row width is too thin for stencil! %D %D",y,s);
306:   ys = 0;
307:   for (i=0; i<(rank % (m*n))/m; i++) { ys += ly[i];}

309:   if (!lz) {
310:     PetscMalloc(p*sizeof(PetscInt), &da->lz);
311:     lz = da->lz;
312:     for (i=0; i<p; i++) {
313:       lz[i] = P/p + ((P % p) > (i % p));
314:     }
315:   }
316:   z  = lz[rank/(m*n)];
317:   if (p > 1 && z < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Plane width is too thin for stencil! %D %D",z,s);
318:   zs = 0;
319:   for (i=0; i<(rank/(m*n)); i++) { zs += lz[i];}
320:   ye = ys + y;
321:   xe = xs + x;
322:   ze = zs + z;

324:   /* determine ghost region */
325:   /* Assume No Periodicity */
326:   if (xs-s > 0) Xs = xs - s; else Xs = 0;
327:   if (ys-s > 0) Ys = ys - s; else Ys = 0;
328:   if (zs-s > 0) Zs = zs - s; else Zs = 0;
329:   if (xe+s <= M) Xe = xe + s; else Xe = M;
330:   if (ye+s <= N) Ye = ye + s; else Ye = N;
331:   if (ze+s <= P) Ze = ze + s; else Ze = P;

333:   /* X Periodic */
334:   if (DAXPeriodic(wrap)){
335:     Xs = xs - s;
336:     Xe = xe + s;
337:   }

339:   /* Y Periodic */
340:   if (DAYPeriodic(wrap)){
341:     Ys = ys - s;
342:     Ye = ye + s;
343:   }

345:   /* Z Periodic */
346:   if (DAZPeriodic(wrap)){
347:     Zs = zs - s;
348:     Ze = ze + s;
349:   }

351:   /* Resize all X parameters to reflect w */
352:   x   *= dof;
353:   xs  *= dof;
354:   xe  *= dof;
355:   Xs  *= dof;
356:   Xe  *= dof;
357:   s_x  = s*dof;
358:   s_y  = s;
359:   s_z  = s;

361:   /* determine starting point of each processor */
362:   nn       = x*y*z;
363:   PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);
364:   MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);
365:   bases[0] = 0;
366:   for (i=1; i<=size; i++) {
367:     bases[i] = ldims[i-1];
368:   }
369:   for (i=1; i<=size; i++) {
370:     bases[i] += bases[i-1];
371:   }

373:   /* allocate the base parallel and sequential vectors */
374:   da->Nlocal = x*y*z;
375:   VecCreateMPIWithArray(comm,da->Nlocal,PETSC_DECIDE,0,&global);
376:   VecSetBlockSize(global,dof);
377:   da->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs);
378:   VecCreateSeqWithArray(MPI_COMM_SELF,da->nlocal,0,&local);
379:   VecSetBlockSize(local,dof);

381:   /* generate appropriate vector scatters */
382:   /* local to global inserts non-ghost point region into global */
383:   VecGetOwnershipRange(global,&start,&end);
384:   ISCreateStride(comm,x*y*z,start,1,&to);

386:   left   = xs - Xs;
387:   bottom = ys - Ys; top = bottom + y;
388:   down   = zs - Zs; up  = down + z;
389:   count  = x*(top-bottom)*(up-down);
390:   PetscMalloc(count*sizeof(PetscInt)/dof,&idx);
391:   count  = 0;
392:   for (i=down; i<up; i++) {
393:     for (j=bottom; j<top; j++) {
394:       for (k=0; k<x; k += dof) {
395:         idx[count++] = (left+j*(Xe-Xs))+i*(Xe-Xs)*(Ye-Ys) + k;
396:       }
397:     }
398:   }

400:   ISCreateBlock(comm,dof,count,idx,&from);
401:   PetscFree(idx);

403:   VecScatterCreate(local,from,global,to,&ltog);
404:   PetscLogObjectParent(da,ltog);
405:   ISDestroy(from);
406:   ISDestroy(to);

408:   /* global to local must include ghost points */
409:   if (stencil_type == DA_STENCIL_BOX) {
410:     ISCreateStride(comm,(Xe-Xs)*(Ye-Ys)*(Ze-Zs),0,1,&to);
411:   } else {
412:     /* This is way ugly! We need to list the funny cross type region */
413:     /* the bottom chunck */
414:     left   = xs - Xs;
415:     bottom = ys - Ys; top = bottom + y;
416:     down   = zs - Zs;   up  = down + z;
417:     count  = down*(top-bottom)*x + (up-down)*(bottom*x  + (top-bottom)*(Xe-Xs) + (Ye-Ys-top)*x) + (Ze-Zs-up)*(top-bottom)*x;
418:     PetscMalloc(count*sizeof(PetscInt)/dof,&idx);
419:     count  = 0;
420:     for (i=0; i<down; i++) {
421:       for (j=bottom; j<top; j++) {
422:         for (k=0; k<x; k += dof) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
423:       }
424:     }
425:     /* the middle piece */
426:     for (i=down; i<up; i++) {
427:       /* front */
428:       for (j=0; j<bottom; j++) {
429:         for (k=0; k<x; k += dof) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
430:       }
431:       /* middle */
432:       for (j=bottom; j<top; j++) {
433:         for (k=0; k<Xe-Xs; k += dof) idx[count++] = j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
434:       }
435:       /* back */
436:       for (j=top; j<Ye-Ys; j++) {
437:         for (k=0; k<x; k += dof) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
438:       }
439:     }
440:     /* the top piece */
441:     for (i=up; i<Ze-Zs; i++) {
442:       for (j=bottom; j<top; j++) {
443:         for (k=0; k<x; k += dof) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
444:       }
445:     }
446:     ISCreateBlock(comm,dof,count,idx,&to);
447:     PetscFree(idx);
448:   }

450:   /* determine who lies on each side of use stored in    n24 n25 n26
451:                                                          n21 n22 n23
452:                                                          n18 n19 n20

454:                                                          n15 n16 n17
455:                                                          n12     n14
456:                                                          n9  n10 n11

458:                                                          n6  n7  n8
459:                                                          n3  n4  n5
460:                                                          n0  n1  n2
461:   */
462: 
463:   /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
464: 
465:   /* Assume Nodes are Internal to the Cube */
466: 
467:   n0  = rank - m*n - m - 1;
468:   n1  = rank - m*n - m;
469:   n2  = rank - m*n - m + 1;
470:   n3  = rank - m*n -1;
471:   n4  = rank - m*n;
472:   n5  = rank - m*n + 1;
473:   n6  = rank - m*n + m - 1;
474:   n7  = rank - m*n + m;
475:   n8  = rank - m*n + m + 1;

477:   n9  = rank - m - 1;
478:   n10 = rank - m;
479:   n11 = rank - m + 1;
480:   n12 = rank - 1;
481:   n14 = rank + 1;
482:   n15 = rank + m - 1;
483:   n16 = rank + m;
484:   n17 = rank + m + 1;

486:   n18 = rank + m*n - m - 1;
487:   n19 = rank + m*n - m;
488:   n20 = rank + m*n - m + 1;
489:   n21 = rank + m*n - 1;
490:   n22 = rank + m*n;
491:   n23 = rank + m*n + 1;
492:   n24 = rank + m*n + m - 1;
493:   n25 = rank + m*n + m;
494:   n26 = rank + m*n + m + 1;

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

498:   if (xs == 0) { /* First assume not corner or edge */
499:     n0  = rank       -1 - (m*n);
500:     n3  = rank + m   -1 - (m*n);
501:     n6  = rank + 2*m -1 - (m*n);
502:     n9  = rank       -1;
503:     n12 = rank + m   -1;
504:     n15 = rank + 2*m -1;
505:     n18 = rank       -1 + (m*n);
506:     n21 = rank + m   -1 + (m*n);
507:     n24 = rank + 2*m -1 + (m*n);
508:    }

510:   if (xe == M*dof) { /* First assume not corner or edge */
511:     n2  = rank -2*m +1 - (m*n);
512:     n5  = rank - m  +1 - (m*n);
513:     n8  = rank      +1 - (m*n);
514:     n11 = rank -2*m +1;
515:     n14 = rank - m  +1;
516:     n17 = rank      +1;
517:     n20 = rank -2*m +1 + (m*n);
518:     n23 = rank - m  +1 + (m*n);
519:     n26 = rank      +1 + (m*n);
520:   }

522:   if (ys==0) { /* First assume not corner or edge */
523:     n0  = rank + m * (n-1) -1 - (m*n);
524:     n1  = rank + m * (n-1)    - (m*n);
525:     n2  = rank + m * (n-1) +1 - (m*n);
526:     n9  = rank + m * (n-1) -1;
527:     n10 = rank + m * (n-1);
528:     n11 = rank + m * (n-1) +1;
529:     n18 = rank + m * (n-1) -1 + (m*n);
530:     n19 = rank + m * (n-1)    + (m*n);
531:     n20 = rank + m * (n-1) +1 + (m*n);
532:   }

534:   if (ye == N) { /* First assume not corner or edge */
535:     n6  = rank - m * (n-1) -1 - (m*n);
536:     n7  = rank - m * (n-1)    - (m*n);
537:     n8  = rank - m * (n-1) +1 - (m*n);
538:     n15 = rank - m * (n-1) -1;
539:     n16 = rank - m * (n-1);
540:     n17 = rank - m * (n-1) +1;
541:     n24 = rank - m * (n-1) -1 + (m*n);
542:     n25 = rank - m * (n-1)    + (m*n);
543:     n26 = rank - m * (n-1) +1 + (m*n);
544:   }
545: 
546:   if (zs == 0) { /* First assume not corner or edge */
547:     n0 = size - (m*n) + rank - m - 1;
548:     n1 = size - (m*n) + rank - m;
549:     n2 = size - (m*n) + rank - m + 1;
550:     n3 = size - (m*n) + rank - 1;
551:     n4 = size - (m*n) + rank;
552:     n5 = size - (m*n) + rank + 1;
553:     n6 = size - (m*n) + rank + m - 1;
554:     n7 = size - (m*n) + rank + m ;
555:     n8 = size - (m*n) + rank + m + 1;
556:   }

558:   if (ze == P) { /* First assume not corner or edge */
559:     n18 = (m*n) - (size-rank) - m - 1;
560:     n19 = (m*n) - (size-rank) - m;
561:     n20 = (m*n) - (size-rank) - m + 1;
562:     n21 = (m*n) - (size-rank) - 1;
563:     n22 = (m*n) - (size-rank);
564:     n23 = (m*n) - (size-rank) + 1;
565:     n24 = (m*n) - (size-rank) + m - 1;
566:     n25 = (m*n) - (size-rank) + m;
567:     n26 = (m*n) - (size-rank) + m + 1;
568:   }

570:   if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */
571:     n0 = size - m*n + rank + m-1 - m;
572:     n3 = size - m*n + rank + m-1;
573:     n6 = size - m*n + rank + m-1 + m;
574:   }
575: 
576:   if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */
577:     n18 = m*n - (size - rank) + m-1 - m;
578:     n21 = m*n - (size - rank) + m-1;
579:     n24 = m*n - (size - rank) + m-1 + m;
580:   }

582:   if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */
583:     n0  = rank + m*n -1 - m*n;
584:     n9  = rank + m*n -1;
585:     n18 = rank + m*n -1 + m*n;
586:   }

588:   if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */
589:     n6  = rank - m*(n-1) + m-1 - m*n;
590:     n15 = rank - m*(n-1) + m-1;
591:     n24 = rank - m*(n-1) + m-1 + m*n;
592:   }

594:   if ((xe==M*dof) && (zs==0)) { /* Assume an edge, not corner */
595:     n2 = size - (m*n-rank) - (m-1) - m;
596:     n5 = size - (m*n-rank) - (m-1);
597:     n8 = size - (m*n-rank) - (m-1) + m;
598:   }

600:   if ((xe==M*dof) && (ze==P)) { /* Assume an edge, not corner */
601:     n20 = m*n - (size - rank) - (m-1) - m;
602:     n23 = m*n - (size - rank) - (m-1);
603:     n26 = m*n - (size - rank) - (m-1) + m;
604:   }

606:   if ((xe==M*dof) && (ys==0)) { /* Assume an edge, not corner */
607:     n2  = rank + m*(n-1) - (m-1) - m*n;
608:     n11 = rank + m*(n-1) - (m-1);
609:     n20 = rank + m*(n-1) - (m-1) + m*n;
610:   }

612:   if ((xe==M*dof) && (ye==N)) { /* Assume an edge, not corner */
613:     n8  = rank - m*n +1 - m*n;
614:     n17 = rank - m*n +1;
615:     n26 = rank - m*n +1 + m*n;
616:   }

618:   if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */
619:     n0 = size - m + rank -1;
620:     n1 = size - m + rank;
621:     n2 = size - m + rank +1;
622:   }

624:   if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */
625:     n18 = m*n - (size - rank) + m*(n-1) -1;
626:     n19 = m*n - (size - rank) + m*(n-1);
627:     n20 = m*n - (size - rank) + m*(n-1) +1;
628:   }

630:   if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */
631:     n6 = size - (m*n-rank) - m * (n-1) -1;
632:     n7 = size - (m*n-rank) - m * (n-1);
633:     n8 = size - (m*n-rank) - m * (n-1) +1;
634:   }

636:   if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */
637:     n24 = rank - (size-m) -1;
638:     n25 = rank - (size-m);
639:     n26 = rank - (size-m) +1;
640:   }

642:   /* Check for Corners */
643:   if ((xs==0)   && (ys==0) && (zs==0)) { n0  = size -1;}
644:   if ((xs==0)   && (ys==0) && (ze==P)) { n18 = m*n-1;}
645:   if ((xs==0)   && (ye==N) && (zs==0)) { n6  = (size-1)-m*(n-1);}
646:   if ((xs==0)   && (ye==N) && (ze==P)) { n24 = m-1;}
647:   if ((xe==M*dof) && (ys==0) && (zs==0)) { n2  = size-m;}
648:   if ((xe==M*dof) && (ys==0) && (ze==P)) { n20 = m*n-m;}
649:   if ((xe==M*dof) && (ye==N) && (zs==0)) { n8  = size-m*n;}
650:   if ((xe==M*dof) && (ye==N) && (ze==P)) { n26 = 0;}

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

654:   /* If not X periodic */
655:   if ((wrap != DA_XPERIODIC)  && (wrap != DA_XYPERIODIC) &&
656:      (wrap != DA_XZPERIODIC) && (wrap != DA_XYZPERIODIC)) {
657:     if (xs==0)   {n0  = n3  = n6  = n9  = n12 = n15 = n18 = n21 = n24 = -2;}
658:     if (xe==M*dof) {n2  = n5  = n8  = n11 = n14 = n17 = n20 = n23 = n26 = -2;}
659:   }

661:   /* If not Y periodic */
662:   if ((wrap != DA_YPERIODIC)  && (wrap != DA_XYPERIODIC) &&
663:       (wrap != DA_YZPERIODIC) && (wrap != DA_XYZPERIODIC)) {
664:     if (ys==0)   {n0  = n1  = n2  = n9  = n10 = n11 = n18 = n19 = n20 = -2;}
665:     if (ye==N)   {n6  = n7  = n8  = n15 = n16 = n17 = n24 = n25 = n26 = -2;}
666:   }

668:   /* If not Z periodic */
669:   if ((wrap != DA_ZPERIODIC)  && (wrap != DA_XZPERIODIC) &&
670:       (wrap != DA_YZPERIODIC) && (wrap != DA_XYZPERIODIC)) {
671:     if (zs==0)   {n0  = n1  = n2  = n3  = n4  = n5  = n6  = n7  = n8  = -2;}
672:     if (ze==P)   {n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;}
673:   }

675:   PetscMalloc(27*sizeof(PetscInt),&da->neighbors);
676:   da->neighbors[0] = n0;
677:   da->neighbors[1] = n1;
678:   da->neighbors[2] = n2;
679:   da->neighbors[3] = n3;
680:   da->neighbors[4] = n4;
681:   da->neighbors[5] = n5;
682:   da->neighbors[6] = n6;
683:   da->neighbors[7] = n7;
684:   da->neighbors[8] = n8;
685:   da->neighbors[9] = n9;
686:   da->neighbors[10] = n10;
687:   da->neighbors[11] = n11;
688:   da->neighbors[12] = n12;
689:   da->neighbors[13] = rank;
690:   da->neighbors[14] = n14;
691:   da->neighbors[15] = n15;
692:   da->neighbors[16] = n16;
693:   da->neighbors[17] = n17;
694:   da->neighbors[18] = n18;
695:   da->neighbors[19] = n19;
696:   da->neighbors[20] = n20;
697:   da->neighbors[21] = n21;
698:   da->neighbors[22] = n22;
699:   da->neighbors[23] = n23;
700:   da->neighbors[24] = n24;
701:   da->neighbors[25] = n25;
702:   da->neighbors[26] = n26;

704:   /* If star stencil then delete the corner neighbors */
705:   if (stencil_type == DA_STENCIL_STAR) {
706:      /* save information about corner neighbors */
707:      sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7;
708:      sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18;
709:      sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25;
710:      sn26 = n26;
711:      n0  = n1  = n2  = n3  = n5  = n6  = n7  = n8  = n9  = n11 =
712:      n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1;
713:   }


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

719:   nn = 0;

721:   /* Bottom Level */
722:   for (k=0; k<s_z; k++) {
723:     for (i=1; i<=s_y; i++) {
724:       if (n0 >= 0) { /* left below */
725:         x_t = lx[n0 % m]*dof;
726:         y_t = ly[(n0 % (m*n))/m];
727:         z_t = lz[n0 / (m*n)];
728:         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;
729:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
730:       }
731:       if (n1 >= 0) { /* directly below */
732:         x_t = x;
733:         y_t = ly[(n1 % (m*n))/m];
734:         z_t = lz[n1 / (m*n)];
735:         s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
736:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
737:       }
738:       if (n2 >= 0) { /* right below */
739:         x_t = lx[n2 % m]*dof;
740:         y_t = ly[(n2 % (m*n))/m];
741:         z_t = lz[n2 / (m*n)];
742:         s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
743:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
744:       }
745:     }

747:     for (i=0; i<y; i++) {
748:       if (n3 >= 0) { /* directly left */
749:         x_t = lx[n3 % m]*dof;
750:         y_t = y;
751:         z_t = lz[n3 / (m*n)];
752:         s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
753:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
754:       }

756:       if (n4 >= 0) { /* middle */
757:         x_t = x;
758:         y_t = y;
759:         z_t = lz[n4 / (m*n)];
760:         s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
761:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
762:       }

764:       if (n5 >= 0) { /* directly right */
765:         x_t = lx[n5 % m]*dof;
766:         y_t = y;
767:         z_t = lz[n5 / (m*n)];
768:         s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
769:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
770:       }
771:     }

773:     for (i=1; i<=s_y; i++) {
774:       if (n6 >= 0) { /* left above */
775:         x_t = lx[n6 % m]*dof;
776:         y_t = ly[(n6 % (m*n))/m];
777:         z_t = lz[n6 / (m*n)];
778:         s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
779:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
780:       }
781:       if (n7 >= 0) { /* directly above */
782:         x_t = x;
783:         y_t = ly[(n7 % (m*n))/m];
784:         z_t = lz[n7 / (m*n)];
785:         s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
786:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
787:       }
788:       if (n8 >= 0) { /* right above */
789:         x_t = lx[n8 % m]*dof;
790:         y_t = ly[(n8 % (m*n))/m];
791:         z_t = lz[n8 / (m*n)];
792:         s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
793:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
794:       }
795:     }
796:   }

798:   /* Middle Level */
799:   for (k=0; k<z; k++) {
800:     for (i=1; i<=s_y; i++) {
801:       if (n9 >= 0) { /* left below */
802:         x_t = lx[n9 % m]*dof;
803:         y_t = ly[(n9 % (m*n))/m];
804:         /* z_t = z; */
805:         s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
806:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
807:       }
808:       if (n10 >= 0) { /* directly below */
809:         x_t = x;
810:         y_t = ly[(n10 % (m*n))/m];
811:         /* z_t = z; */
812:         s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
813:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
814:       }
815:       if (n11 >= 0) { /* right below */
816:         x_t = lx[n11 % m]*dof;
817:         y_t = ly[(n11 % (m*n))/m];
818:         /* z_t = z; */
819:         s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
820:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
821:       }
822:     }

824:     for (i=0; i<y; i++) {
825:       if (n12 >= 0) { /* directly left */
826:         x_t = lx[n12 % m]*dof;
827:         y_t = y;
828:         /* z_t = z; */
829:         s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
830:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
831:       }

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

837:       if (n14 >= 0) { /* directly right */
838:         x_t = lx[n14 % m]*dof;
839:         y_t = y;
840:         /* z_t = z; */
841:         s_t = bases[n14] + i*x_t + k*x_t*y_t;
842:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
843:       }
844:     }

846:     for (i=1; i<=s_y; i++) {
847:       if (n15 >= 0) { /* left above */
848:         x_t = lx[n15 % m]*dof;
849:         y_t = ly[(n15 % (m*n))/m];
850:         /* z_t = z; */
851:         s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
852:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
853:       }
854:       if (n16 >= 0) { /* directly above */
855:         x_t = x;
856:         y_t = ly[(n16 % (m*n))/m];
857:         /* z_t = z; */
858:         s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
859:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
860:       }
861:       if (n17 >= 0) { /* right above */
862:         x_t = lx[n17 % m]*dof;
863:         y_t = ly[(n17 % (m*n))/m];
864:         /* z_t = z; */
865:         s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
866:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
867:       }
868:     }
869:   }
870: 
871:   /* Upper Level */
872:   for (k=0; k<s_z; k++) {
873:     for (i=1; i<=s_y; i++) {
874:       if (n18 >= 0) { /* left below */
875:         x_t = lx[n18 % m]*dof;
876:         y_t = ly[(n18 % (m*n))/m];
877:         /* z_t = lz[n18 / (m*n)]; */
878:         s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
879:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
880:       }
881:       if (n19 >= 0) { /* directly below */
882:         x_t = x;
883:         y_t = ly[(n19 % (m*n))/m];
884:         /* z_t = lz[n19 / (m*n)]; */
885:         s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
886:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
887:       }
888:       if (n20 >= 0) { /* right below */
889:         x_t = lx[n20 % m]*dof;
890:         y_t = ly[(n20 % (m*n))/m];
891:         /* z_t = lz[n20 / (m*n)]; */
892:         s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
893:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
894:       }
895:     }

897:     for (i=0; i<y; i++) {
898:       if (n21 >= 0) { /* directly left */
899:         x_t = lx[n21 % m]*dof;
900:         y_t = y;
901:         /* z_t = lz[n21 / (m*n)]; */
902:         s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
903:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
904:       }

906:       if (n22 >= 0) { /* middle */
907:         x_t = x;
908:         y_t = y;
909:         /* z_t = lz[n22 / (m*n)]; */
910:         s_t = bases[n22] + i*x_t + k*x_t*y_t;
911:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
912:       }

914:       if (n23 >= 0) { /* directly right */
915:         x_t = lx[n23 % m]*dof;
916:         y_t = y;
917:         /* z_t = lz[n23 / (m*n)]; */
918:         s_t = bases[n23] + i*x_t + k*x_t*y_t;
919:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
920:       }
921:     }

923:     for (i=1; i<=s_y; i++) {
924:       if (n24 >= 0) { /* left above */
925:         x_t = lx[n24 % m]*dof;
926:         y_t = ly[(n24 % (m*n))/m];
927:         /* z_t = lz[n24 / (m*n)]; */
928:         s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
929:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
930:       }
931:       if (n25 >= 0) { /* directly above */
932:         x_t = x;
933:         y_t = ly[(n25 % (m*n))/m];
934:         /* z_t = lz[n25 / (m*n)]; */
935:         s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
936:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
937:       }
938:       if (n26 >= 0) { /* right above */
939:         x_t = lx[n26 % m]*dof;
940:         y_t = ly[(n26 % (m*n))/m];
941:         /* z_t = lz[n26 / (m*n)]; */
942:         s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
943:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
944:       }
945:     }
946:   }
947:   base = bases[rank];
948:   {
949:     PetscInt nnn = nn/dof,*iidx;
950:     PetscMalloc(nnn*sizeof(PetscInt),&iidx);
951:     for (i=0; i<nnn; i++) {
952:       iidx[i] = idx[dof*i];
953:     }
954:     ISCreateBlock(comm,dof,nnn,iidx,&from);
955:     PetscFree(iidx);
956:   }
957:   VecScatterCreate(global,from,local,to,&gtol);
958:   PetscLogObjectParent(da,gtol);
959:   ISDestroy(from);
960:   ISDestroy(to);

962:   da->m  = m;  da->n  = n;  da->p  = p;
963:   da->xs = xs; da->xe = xe; da->ys = ys; da->ye = ye; da->zs = zs; da->ze = ze;
964:   da->Xs = Xs; da->Xe = Xe; da->Ys = Ys; da->Ye = Ye; da->Zs = Zs; da->Ze = Ze;

966:   VecDestroy(local);
967:   VecDestroy(global);

969:   if (stencil_type == DA_STENCIL_STAR) {
970:     /*
971:         Recompute the local to global mappings, this time keeping the 
972:       information about the cross corner processor numbers.
973:     */
974:     n0  = sn0;  n1  = sn1;  n2  = sn2;  n3  = sn3;  n5  = sn5;  n6  = sn6; n7 = sn7;
975:     n8  = sn8;  n9  = sn9;  n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18;
976:     n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25;
977:     n26 = sn26;

979:     nn = 0;

981:     /* Bottom Level */
982:     for (k=0; k<s_z; k++) {
983:       for (i=1; i<=s_y; i++) {
984:         if (n0 >= 0) { /* left below */
985:           x_t = lx[n0 % m]*dof;
986:           y_t = ly[(n0 % (m*n))/m];
987:           z_t = lz[n0 / (m*n)];
988:           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;
989:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
990:         }
991:         if (n1 >= 0) { /* directly below */
992:           x_t = x;
993:           y_t = ly[(n1 % (m*n))/m];
994:           z_t = lz[n1 / (m*n)];
995:           s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
996:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
997:         }
998:         if (n2 >= 0) { /* right below */
999:           x_t = lx[n2 % m]*dof;
1000:           y_t = ly[(n2 % (m*n))/m];
1001:           z_t = lz[n2 / (m*n)];
1002:           s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
1003:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1004:         }
1005:       }

1007:       for (i=0; i<y; i++) {
1008:         if (n3 >= 0) { /* directly left */
1009:           x_t = lx[n3 % m]*dof;
1010:           y_t = y;
1011:           z_t = lz[n3 / (m*n)];
1012:           s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1013:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1014:         }

1016:         if (n4 >= 0) { /* middle */
1017:           x_t = x;
1018:           y_t = y;
1019:           z_t = lz[n4 / (m*n)];
1020:           s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1021:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1022:         }

1024:         if (n5 >= 0) { /* directly right */
1025:           x_t = lx[n5 % m]*dof;
1026:           y_t = y;
1027:           z_t = lz[n5 / (m*n)];
1028:           s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1029:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1030:         }
1031:       }

1033:       for (i=1; i<=s_y; i++) {
1034:         if (n6 >= 0) { /* left above */
1035:           x_t = lx[n6 % m]*dof;
1036:           y_t = ly[(n6 % (m*n))/m];
1037:           z_t = lz[n6 / (m*n)];
1038:           s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1039:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1040:         }
1041:         if (n7 >= 0) { /* directly above */
1042:           x_t = x;
1043:           y_t = ly[(n7 % (m*n))/m];
1044:           z_t = lz[n7 / (m*n)];
1045:           s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1046:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1047:         }
1048:         if (n8 >= 0) { /* right above */
1049:           x_t = lx[n8 % m]*dof;
1050:           y_t = ly[(n8 % (m*n))/m];
1051:           z_t = lz[n8 / (m*n)];
1052:           s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1053:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1054:         }
1055:       }
1056:     }

1058:     /* Middle Level */
1059:     for (k=0; k<z; k++) {
1060:       for (i=1; i<=s_y; i++) {
1061:         if (n9 >= 0) { /* left below */
1062:           x_t = lx[n9 % m]*dof;
1063:           y_t = ly[(n9 % (m*n))/m];
1064:           /* z_t = z; */
1065:           s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
1066:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1067:         }
1068:         if (n10 >= 0) { /* directly below */
1069:           x_t = x;
1070:           y_t = ly[(n10 % (m*n))/m];
1071:           /* z_t = z; */
1072:           s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1073:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1074:         }
1075:         if (n11 >= 0) { /* right below */
1076:           x_t = lx[n11 % m]*dof;
1077:           y_t = ly[(n11 % (m*n))/m];
1078:           /* z_t = z; */
1079:           s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1080:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1081:         }
1082:       }

1084:       for (i=0; i<y; i++) {
1085:         if (n12 >= 0) { /* directly left */
1086:           x_t = lx[n12 % m]*dof;
1087:           y_t = y;
1088:           /* z_t = z; */
1089:           s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
1090:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1091:         }

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

1097:         if (n14 >= 0) { /* directly right */
1098:           x_t = lx[n14 % m]*dof;
1099:           y_t = y;
1100:           /* z_t = z; */
1101:           s_t = bases[n14] + i*x_t + k*x_t*y_t;
1102:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1103:         }
1104:       }

1106:       for (i=1; i<=s_y; i++) {
1107:         if (n15 >= 0) { /* left above */
1108:           x_t = lx[n15 % m]*dof;
1109:           y_t = ly[(n15 % (m*n))/m];
1110:           /* z_t = z; */
1111:           s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
1112:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1113:         }
1114:         if (n16 >= 0) { /* directly above */
1115:           x_t = x;
1116:           y_t = ly[(n16 % (m*n))/m];
1117:           /* z_t = z; */
1118:           s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
1119:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1120:         }
1121:         if (n17 >= 0) { /* right above */
1122:           x_t = lx[n17 % m]*dof;
1123:           y_t = ly[(n17 % (m*n))/m];
1124:           /* z_t = z; */
1125:           s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
1126:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1127:         }
1128:       }
1129:     }
1130: 
1131:     /* Upper Level */
1132:     for (k=0; k<s_z; k++) {
1133:       for (i=1; i<=s_y; i++) {
1134:         if (n18 >= 0) { /* left below */
1135:           x_t = lx[n18 % m]*dof;
1136:           y_t = ly[(n18 % (m*n))/m];
1137:           /* z_t = lz[n18 / (m*n)]; */
1138:           s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
1139:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1140:         }
1141:         if (n19 >= 0) { /* directly below */
1142:           x_t = x;
1143:           y_t = ly[(n19 % (m*n))/m];
1144:           /* z_t = lz[n19 / (m*n)]; */
1145:           s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1146:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1147:         }
1148:         if (n20 >= 0) { /* right below */
1149:           x_t = lx[n20 % m]*dof;
1150:           y_t = ly[(n20 % (m*n))/m];
1151:           /* z_t = lz[n20 / (m*n)]; */
1152:           s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1153:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1154:         }
1155:       }

1157:       for (i=0; i<y; i++) {
1158:         if (n21 >= 0) { /* directly left */
1159:           x_t = lx[n21 % m]*dof;
1160:           y_t = y;
1161:           /* z_t = lz[n21 / (m*n)]; */
1162:           s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
1163:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1164:         }

1166:         if (n22 >= 0) { /* middle */
1167:           x_t = x;
1168:           y_t = y;
1169:           /* z_t = lz[n22 / (m*n)]; */
1170:           s_t = bases[n22] + i*x_t + k*x_t*y_t;
1171:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1172:         }

1174:         if (n23 >= 0) { /* directly right */
1175:           x_t = lx[n23 % m]*dof;
1176:           y_t = y;
1177:           /* z_t = lz[n23 / (m*n)]; */
1178:           s_t = bases[n23] + i*x_t + k*x_t*y_t;
1179:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1180:         }
1181:       }

1183:       for (i=1; i<=s_y; i++) {
1184:         if (n24 >= 0) { /* left above */
1185:           x_t = lx[n24 % m]*dof;
1186:           y_t = ly[(n24 % (m*n))/m];
1187:           /* z_t = lz[n24 / (m*n)]; */
1188:           s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
1189:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1190:         }
1191:         if (n25 >= 0) { /* directly above */
1192:           x_t = x;
1193:           y_t = ly[(n25 % (m*n))/m];
1194:           /* z_t = lz[n25 / (m*n)]; */
1195:           s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
1196:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1197:         }
1198:         if (n26 >= 0) { /* right above */
1199:           x_t = lx[n26 % m]*dof;
1200:           y_t = ly[(n26 % (m*n))/m];
1201:           /* z_t = lz[n26 / (m*n)]; */
1202:           s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
1203:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1204:         }
1205:       }
1206:     }
1207:   }
1208:   /* redo idx to include "missing" ghost points */
1209:   /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
1210: 
1211:   /* Assume Nodes are Internal to the Cube */
1212: 
1213:   n0  = rank - m*n - m - 1;
1214:   n1  = rank - m*n - m;
1215:   n2  = rank - m*n - m + 1;
1216:   n3  = rank - m*n -1;
1217:   n4  = rank - m*n;
1218:   n5  = rank - m*n + 1;
1219:   n6  = rank - m*n + m - 1;
1220:   n7  = rank - m*n + m;
1221:   n8  = rank - m*n + m + 1;

1223:   n9  = rank - m - 1;
1224:   n10 = rank - m;
1225:   n11 = rank - m + 1;
1226:   n12 = rank - 1;
1227:   n14 = rank + 1;
1228:   n15 = rank + m - 1;
1229:   n16 = rank + m;
1230:   n17 = rank + m + 1;

1232:   n18 = rank + m*n - m - 1;
1233:   n19 = rank + m*n - m;
1234:   n20 = rank + m*n - m + 1;
1235:   n21 = rank + m*n - 1;
1236:   n22 = rank + m*n;
1237:   n23 = rank + m*n + 1;
1238:   n24 = rank + m*n + m - 1;
1239:   n25 = rank + m*n + m;
1240:   n26 = rank + m*n + m + 1;

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

1244:   if (xs == 0) { /* First assume not corner or edge */
1245:     n0  = rank       -1 - (m*n);
1246:     n3  = rank + m   -1 - (m*n);
1247:     n6  = rank + 2*m -1 - (m*n);
1248:     n9  = rank       -1;
1249:     n12 = rank + m   -1;
1250:     n15 = rank + 2*m -1;
1251:     n18 = rank       -1 + (m*n);
1252:     n21 = rank + m   -1 + (m*n);
1253:     n24 = rank + 2*m -1 + (m*n);
1254:    }

1256:   if (xe == M*dof) { /* First assume not corner or edge */
1257:     n2  = rank -2*m +1 - (m*n);
1258:     n5  = rank - m  +1 - (m*n);
1259:     n8  = rank      +1 - (m*n);
1260:     n11 = rank -2*m +1;
1261:     n14 = rank - m  +1;
1262:     n17 = rank      +1;
1263:     n20 = rank -2*m +1 + (m*n);
1264:     n23 = rank - m  +1 + (m*n);
1265:     n26 = rank      +1 + (m*n);
1266:   }

1268:   if (ys==0) { /* First assume not corner or edge */
1269:     n0  = rank + m * (n-1) -1 - (m*n);
1270:     n1  = rank + m * (n-1)    - (m*n);
1271:     n2  = rank + m * (n-1) +1 - (m*n);
1272:     n9  = rank + m * (n-1) -1;
1273:     n10 = rank + m * (n-1);
1274:     n11 = rank + m * (n-1) +1;
1275:     n18 = rank + m * (n-1) -1 + (m*n);
1276:     n19 = rank + m * (n-1)    + (m*n);
1277:     n20 = rank + m * (n-1) +1 + (m*n);
1278:   }

1280:   if (ye == N) { /* First assume not corner or edge */
1281:     n6  = rank - m * (n-1) -1 - (m*n);
1282:     n7  = rank - m * (n-1)    - (m*n);
1283:     n8  = rank - m * (n-1) +1 - (m*n);
1284:     n15 = rank - m * (n-1) -1;
1285:     n16 = rank - m * (n-1);
1286:     n17 = rank - m * (n-1) +1;
1287:     n24 = rank - m * (n-1) -1 + (m*n);
1288:     n25 = rank - m * (n-1)    + (m*n);
1289:     n26 = rank - m * (n-1) +1 + (m*n);
1290:   }
1291: 
1292:   if (zs == 0) { /* First assume not corner or edge */
1293:     n0 = size - (m*n) + rank - m - 1;
1294:     n1 = size - (m*n) + rank - m;
1295:     n2 = size - (m*n) + rank - m + 1;
1296:     n3 = size - (m*n) + rank - 1;
1297:     n4 = size - (m*n) + rank;
1298:     n5 = size - (m*n) + rank + 1;
1299:     n6 = size - (m*n) + rank + m - 1;
1300:     n7 = size - (m*n) + rank + m ;
1301:     n8 = size - (m*n) + rank + m + 1;
1302:   }

1304:   if (ze == P) { /* First assume not corner or edge */
1305:     n18 = (m*n) - (size-rank) - m - 1;
1306:     n19 = (m*n) - (size-rank) - m;
1307:     n20 = (m*n) - (size-rank) - m + 1;
1308:     n21 = (m*n) - (size-rank) - 1;
1309:     n22 = (m*n) - (size-rank);
1310:     n23 = (m*n) - (size-rank) + 1;
1311:     n24 = (m*n) - (size-rank) + m - 1;
1312:     n25 = (m*n) - (size-rank) + m;
1313:     n26 = (m*n) - (size-rank) + m + 1;
1314:   }

1316:   if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */
1317:     n0 = size - m*n + rank + m-1 - m;
1318:     n3 = size - m*n + rank + m-1;
1319:     n6 = size - m*n + rank + m-1 + m;
1320:   }
1321: 
1322:   if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */
1323:     n18 = m*n - (size - rank) + m-1 - m;
1324:     n21 = m*n - (size - rank) + m-1;
1325:     n24 = m*n - (size - rank) + m-1 + m;
1326:   }

1328:   if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */
1329:     n0  = rank + m*n -1 - m*n;
1330:     n9  = rank + m*n -1;
1331:     n18 = rank + m*n -1 + m*n;
1332:   }

1334:   if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */
1335:     n6  = rank - m*(n-1) + m-1 - m*n;
1336:     n15 = rank - m*(n-1) + m-1;
1337:     n24 = rank - m*(n-1) + m-1 + m*n;
1338:   }

1340:   if ((xe==M*dof) && (zs==0)) { /* Assume an edge, not corner */
1341:     n2 = size - (m*n-rank) - (m-1) - m;
1342:     n5 = size - (m*n-rank) - (m-1);
1343:     n8 = size - (m*n-rank) - (m-1) + m;
1344:   }

1346:   if ((xe==M*dof) && (ze==P)) { /* Assume an edge, not corner */
1347:     n20 = m*n - (size - rank) - (m-1) - m;
1348:     n23 = m*n - (size - rank) - (m-1);
1349:     n26 = m*n - (size - rank) - (m-1) + m;
1350:   }

1352:   if ((xe==M*dof) && (ys==0)) { /* Assume an edge, not corner */
1353:     n2  = rank + m*(n-1) - (m-1) - m*n;
1354:     n11 = rank + m*(n-1) - (m-1);
1355:     n20 = rank + m*(n-1) - (m-1) + m*n;
1356:   }

1358:   if ((xe==M*dof) && (ye==N)) { /* Assume an edge, not corner */
1359:     n8  = rank - m*n +1 - m*n;
1360:     n17 = rank - m*n +1;
1361:     n26 = rank - m*n +1 + m*n;
1362:   }

1364:   if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */
1365:     n0 = size - m + rank -1;
1366:     n1 = size - m + rank;
1367:     n2 = size - m + rank +1;
1368:   }

1370:   if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */
1371:     n18 = m*n - (size - rank) + m*(n-1) -1;
1372:     n19 = m*n - (size - rank) + m*(n-1);
1373:     n20 = m*n - (size - rank) + m*(n-1) +1;
1374:   }

1376:   if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */
1377:     n6 = size - (m*n-rank) - m * (n-1) -1;
1378:     n7 = size - (m*n-rank) - m * (n-1);
1379:     n8 = size - (m*n-rank) - m * (n-1) +1;
1380:   }

1382:   if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */
1383:     n24 = rank - (size-m) -1;
1384:     n25 = rank - (size-m);
1385:     n26 = rank - (size-m) +1;
1386:   }

1388:   /* Check for Corners */
1389:   if ((xs==0)   && (ys==0) && (zs==0)) { n0  = size -1;}
1390:   if ((xs==0)   && (ys==0) && (ze==P)) { n18 = m*n-1;}
1391:   if ((xs==0)   && (ye==N) && (zs==0)) { n6  = (size-1)-m*(n-1);}
1392:   if ((xs==0)   && (ye==N) && (ze==P)) { n24 = m-1;}
1393:   if ((xe==M*dof) && (ys==0) && (zs==0)) { n2  = size-m;}
1394:   if ((xe==M*dof) && (ys==0) && (ze==P)) { n20 = m*n-m;}
1395:   if ((xe==M*dof) && (ye==N) && (zs==0)) { n8  = size-m*n;}
1396:   if ((xe==M*dof) && (ye==N) && (ze==P)) { n26 = 0;}

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

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

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

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

1418:   nn = 0;

1420:   /* Bottom Level */
1421:   for (k=0; k<s_z; k++) {
1422:     for (i=1; i<=s_y; i++) {
1423:       if (n0 >= 0) { /* left below */
1424:         x_t = lx[n0 % m]*dof;
1425:         y_t = ly[(n0 % (m*n))/m];
1426:         z_t = lz[n0 / (m*n)];
1427:         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;
1428:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1429:       }
1430:       if (n1 >= 0) { /* directly below */
1431:         x_t = x;
1432:         y_t = ly[(n1 % (m*n))/m];
1433:         z_t = lz[n1 / (m*n)];
1434:         s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
1435:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1436:       }
1437:       if (n2 >= 0) { /* right below */
1438:         x_t = lx[n2 % m]*dof;
1439:         y_t = ly[(n2 % (m*n))/m];
1440:         z_t = lz[n2 / (m*n)];
1441:         s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
1442:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1443:       }
1444:     }

1446:     for (i=0; i<y; i++) {
1447:       if (n3 >= 0) { /* directly left */
1448:         x_t = lx[n3 % m]*dof;
1449:         y_t = y;
1450:         z_t = lz[n3 / (m*n)];
1451:         s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1452:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1453:       }

1455:       if (n4 >= 0) { /* middle */
1456:         x_t = x;
1457:         y_t = y;
1458:         z_t = lz[n4 / (m*n)];
1459:         s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1460:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1461:       }

1463:       if (n5 >= 0) { /* directly right */
1464:         x_t = lx[n5 % m]*dof;
1465:         y_t = y;
1466:         z_t = lz[n5 / (m*n)];
1467:         s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1468:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1469:       }
1470:     }

1472:     for (i=1; i<=s_y; i++) {
1473:       if (n6 >= 0) { /* left above */
1474:         x_t = lx[n6 % m]*dof;
1475:         y_t = ly[(n6 % (m*n))/m];
1476:         z_t = lz[n6 / (m*n)];
1477:         s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1478:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1479:       }
1480:       if (n7 >= 0) { /* directly above */
1481:         x_t = x;
1482:         y_t = ly[(n7 % (m*n))/m];
1483:         z_t = lz[n7 / (m*n)];
1484:         s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1485:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1486:       }
1487:       if (n8 >= 0) { /* right above */
1488:         x_t = lx[n8 % m]*dof;
1489:         y_t = ly[(n8 % (m*n))/m];
1490:         z_t = lz[n8 / (m*n)];
1491:         s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1492:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1493:       }
1494:     }
1495:   }

1497:   /* Middle Level */
1498:   for (k=0; k<z; k++) {
1499:     for (i=1; i<=s_y; i++) {
1500:       if (n9 >= 0) { /* left below */
1501:         x_t = lx[n9 % m]*dof;
1502:         y_t = ly[(n9 % (m*n))/m];
1503:         /* z_t = z; */
1504:         s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
1505:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1506:       }
1507:       if (n10 >= 0) { /* directly below */
1508:         x_t = x;
1509:         y_t = ly[(n10 % (m*n))/m];
1510:         /* z_t = z; */
1511:         s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1512:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1513:       }
1514:       if (n11 >= 0) { /* right below */
1515:         x_t = lx[n11 % m]*dof;
1516:         y_t = ly[(n11 % (m*n))/m];
1517:         /* z_t = z; */
1518:         s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1519:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1520:       }
1521:     }

1523:     for (i=0; i<y; i++) {
1524:       if (n12 >= 0) { /* directly left */
1525:         x_t = lx[n12 % m]*dof;
1526:         y_t = y;
1527:         /* z_t = z; */
1528:         s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
1529:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1530:       }

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

1536:       if (n14 >= 0) { /* directly right */
1537:         x_t = lx[n14 % m]*dof;
1538:         y_t = y;
1539:         /* z_t = z; */
1540:         s_t = bases[n14] + i*x_t + k*x_t*y_t;
1541:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1542:       }
1543:     }

1545:     for (i=1; i<=s_y; i++) {
1546:       if (n15 >= 0) { /* left above */
1547:         x_t = lx[n15 % m]*dof;
1548:         y_t = ly[(n15 % (m*n))/m];
1549:         /* z_t = z; */
1550:         s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
1551:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1552:       }
1553:       if (n16 >= 0) { /* directly above */
1554:         x_t = x;
1555:         y_t = ly[(n16 % (m*n))/m];
1556:         /* z_t = z; */
1557:         s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
1558:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1559:       }
1560:       if (n17 >= 0) { /* right above */
1561:         x_t = lx[n17 % m]*dof;
1562:         y_t = ly[(n17 % (m*n))/m];
1563:         /* z_t = z; */
1564:         s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
1565:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1566:       }
1567:     }
1568:   }
1569: 
1570:   /* Upper Level */
1571:   for (k=0; k<s_z; k++) {
1572:     for (i=1; i<=s_y; i++) {
1573:       if (n18 >= 0) { /* left below */
1574:         x_t = lx[n18 % m]*dof;
1575:         y_t = ly[(n18 % (m*n))/m];
1576:         /* z_t = lz[n18 / (m*n)]; */
1577:         s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
1578:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1579:       }
1580:       if (n19 >= 0) { /* directly below */
1581:         x_t = x;
1582:         y_t = ly[(n19 % (m*n))/m];
1583:         /* z_t = lz[n19 / (m*n)]; */
1584:         s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1585:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1586:       }
1587:       if (n20 >= 0) { /* right belodof */
1588:         x_t = lx[n20 % m]*dof;
1589:         y_t = ly[(n20 % (m*n))/m];
1590:         /* z_t = lz[n20 / (m*n)]; */
1591:         s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1592:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1593:       }
1594:     }

1596:     for (i=0; i<y; i++) {
1597:       if (n21 >= 0) { /* directly left */
1598:         x_t = lx[n21 % m]*dof;
1599:         y_t = y;
1600:         /* z_t = lz[n21 / (m*n)]; */
1601:         s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
1602:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1603:       }

1605:       if (n22 >= 0) { /* middle */
1606:         x_t = x;
1607:         y_t = y;
1608:         /* z_t = lz[n22 / (m*n)]; */
1609:         s_t = bases[n22] + i*x_t + k*x_t*y_t;
1610:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1611:       }

1613:       if (n23 >= 0) { /* directly right */
1614:         x_t = lx[n23 % m]*dof;
1615:         y_t = y;
1616:         /* z_t = lz[n23 / (m*n)]; */
1617:         s_t = bases[n23] + i*x_t + k*x_t*y_t;
1618:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1619:       }
1620:     }

1622:     for (i=1; i<=s_y; i++) {
1623:       if (n24 >= 0) { /* left above */
1624:         x_t = lx[n24 % m]*dof;
1625:         y_t = ly[(n24 % (m*n))/m];
1626:         /* z_t = lz[n24 / (m*n)]; */
1627:         s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
1628:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1629:       }
1630:       if (n25 >= 0) { /* directly above */
1631:         x_t = x;
1632:         y_t = ly[(n25 % (m*n))/m];
1633:         /* z_t = lz[n25 / (m*n)]; */
1634:         s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
1635:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1636:       }
1637:       if (n26 >= 0) { /* right above */
1638:         x_t = lx[n26 % m]*dof;
1639:         y_t = ly[(n26 % (m*n))/m];
1640:         /* z_t = lz[n26 / (m*n)]; */
1641:         s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
1642:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1643:       }
1644:     }
1645:   }
1646:   PetscFree2(bases,ldims);
1647:   da->gtol      = gtol;
1648:   da->ltog      = ltog;
1649:   da->idx       = idx;
1650:   da->Nl        = nn;
1651:   da->base      = base;
1652:   da->ops->view = DAView_3d;

1654:   /* 
1655:      Set the local to global ordering in the global vector, this allows use
1656:      of VecSetValuesLocal().
1657:   */
1658:   ISLocalToGlobalMappingCreateNC(comm,nn,idx,&da->ltogmap);
1659:   ISLocalToGlobalMappingBlock(da->ltogmap,da->w,&da->ltogmapb);
1660:   PetscLogObjectParent(da,da->ltogmap);

1662:   da->ltol = PETSC_NULL;
1663:   da->ao   = PETSC_NULL;
1664:   return(0);
1665: }

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

1674:    Collective on MPI_Comm

1676:    Input Parameters:
1677: +  comm - MPI communicator
1678: .  wrap - type of periodicity the array should have, if any.  Use one
1679:           of DA_NONPERIODIC, DA_XPERIODIC, DA_YPERIODIC, DA_ZPERIODIC, DA_XYPERIODIC, DA_XZPERIODIC, DA_YZPERIODIC, DA_XYZPERIODIC, or DA_XYZGHOSTED.
1680: .  stencil_type - Type of stencil (DA_STENCIL_STAR or DA_STENCIL_BOX)
1681: .  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 
1682:             from the command line with -da_grid_x <M> -da_grid_y <N> -da_grid_z <P>)
1683: .  m,n,p - corresponding number of processors in each dimension 
1684:            (or PETSC_DECIDE to have calculated)
1685: .  dof - number of degrees of freedom per node
1686: .  lx, ly, lz - arrays containing the number of nodes in each cell along
1687:           the x, y, and z coordinates, or PETSC_NULL. If non-null, these
1688:           must be of length as m,n,p and the corresponding
1689:           m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of
1690:           the ly[] must N, sum of the lz[] must be P
1691: -  s - stencil width

1693:    Output Parameter:
1694: .  da - the resulting distributed array object

1696:    Options Database Key:
1697: +  -da_view - Calls DAView() at the conclusion of DACreate3d()
1698: .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
1699: .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
1700: .  -da_grid_z <nz> - number of grid points in z direction, if P < 0
1701: .  -da_processors_x <MX> - number of processors in x direction
1702: .  -da_processors_y <MY> - number of processors in y direction
1703: .  -da_processors_z <MZ> - number of processors in z direction
1704: .  -da_refine_x - refinement ratio in x direction
1705: .  -da_refine_y - refinement ratio in y direction
1706: -  -da_refine_y - refinement ratio in z direction

1708:    Level: beginner

1710:    Notes:
1711:    The stencil type DA_STENCIL_STAR with width 1 corresponds to the 
1712:    standard 7-pt stencil, while DA_STENCIL_BOX with width 1 denotes
1713:    the standard 27-pt stencil.

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

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

1721: .seealso: DADestroy(), DAView(), DACreate1d(), DACreate2d(), DAGlobalToLocalBegin(), DAGetRefinementFactor(),
1722:           DAGlobalToLocalEnd(), DALocalToGlobal(), DALocalToLocalBegin(), DALocalToLocalEnd(), DASetRefinementFactor(),
1723:           DAGetInfo(), DACreateGlobalVector(), DACreateLocalVector(), DACreateNaturalVector(), DALoad(), DAView(), DAGetOwnershipRanges()

1725: @*/
1726: PetscErrorCode  DACreate3d(MPI_Comm comm,DAPeriodicType wrap,DAStencilType stencil_type,PetscInt M,
1727:                PetscInt N,PetscInt P,PetscInt m,PetscInt n,PetscInt p,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],const PetscInt lz[],DA *da)
1728: {

1732:   DACreate(comm, da);
1733:   DASetDim(*da, 3);
1734:   DASetSizes(*da, M, N, P);
1735:   DASetNumProcs(*da, m, n, p);
1736:   DASetPeriodicity(*da, wrap);
1737:   DASetDof(*da, dof);
1738:   DASetStencilType(*da, stencil_type);
1739:   DASetStencilWidth(*da, s);
1740:   DASetVertexDivision(*da, lx, ly, lz);
1741:   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
1742:   DASetFromOptions(*da);
1743:   DASetType(*da, DA3D);
1744:   return(0);
1745: }