Actual source code: da3.c

  1: /*$Id: da3.c,v 1.136 2001/09/07 20:12:17 bsmith Exp $*/

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

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

 10: #if defined (PETSC_HAVE_AMS)
 11: EXTERN_C_BEGIN
 12: EXTERN int AMSSetFieldBlock_DA(AMS_Memory,char *,Vec);
 13: EXTERN_C_END
 14: #endif

 16: int DAView_3d(DA da,PetscViewer viewer)
 17: {
 18:   int        rank,ierr;
 19:   PetscTruth isascii,isdraw,isbinary;

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

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

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

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

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

 83:         PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);
 84:         PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);
 85:         PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);
 86:         PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);

 88:         xmin = da->xs/da->w;
 89:         xmax =(da->xe-1)/da->w;

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

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

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

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

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

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

151: EXTERN int DAPublish_Petsc(PetscObject);

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

157:    Collective on MPI_Comm

159:    Input Parameters:
160: +  comm - MPI communicator
161: .  wrap - type of periodicity the array should have, if any.  Use one
162:           of DA_NONPERIODIC, DA_XPERIODIC, DA_YPERIODIC, DA_XYPERIODIC, DA_XYZPERIODIC, DA_XZPERIODIC, or DA_YZPERIODIC.
163: .  stencil_type - Type of stencil (DA_STENCIL_STAR or DA_STENCIL_BOX)
164: .  M,N,P - global dimension in each direction of the array
165: .  m,n,p - corresponding number of processors in each dimension 
166:            (or PETSC_DECIDE to have calculated)
167: .  dof - number of degrees of freedom per node
168: .  lx, ly, lz - arrays containing the number of nodes in each cell along
169:           the x, y, and z coordinates, or PETSC_NULL. If non-null, these
170:           must be of length as m,n,p and the corresponding
171:           m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of
172:           the ly[] must n, sum of the lz[] must be P
173: -  s - stencil width

175:    Output Parameter:
176: .  inra - the resulting distributed array object

178:    Options Database Key:
179: +  -da_view - Calls DAView() at the conclusion of DACreate3d()
180: .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
181: .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
182: -  -da_grid_z <nz> - number of grid points in z direction, if P < 0

184:    Level: beginner

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

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

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

197: .seealso: DADestroy(), DAView(), DACreate1d(), DACreate2d(), DAGlobalToLocalBegin(),
198:           DAGlobalToLocalEnd(), DALocalToGlobal(), DALocalToLocalBegin(), DALocalToLocalEnd(),
199:           DAGetInfo(), DACreateGlobalVector(), DACreateLocalVector(), DACreateNaturalVector(), DALoad(), DAView()

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

465:   /* allocate the base parallel and sequential vectors */
466:   VecCreateMPI(comm,x*y*z,PETSC_DECIDE,&global);
467:   VecSetBlockSize(global,dof);
468:   VecCreateSeq(MPI_COMM_SELF,(Xe-Xs)*(Ye-Ys)*(Ze-Zs),&local);
469:   VecSetBlockSize(local,dof);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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


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

783:   nn = 0;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1026:   PetscLogObjectParent(da,global);
1027:   PetscLogObjectParent(da,local);

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

1039:     nn = 0;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1478:   nn = 0;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1728:   da->ltol = PETSC_NULL;
1729:   da->ao   = PETSC_NULL;

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


1748:   PetscOptionsHasName(PETSC_NULL,"-da_view",&flg1);
1749:   if (flg1) {DAView(da,PETSC_VIEWER_STDOUT_(da->comm));}
1750:   PetscOptionsHasName(PETSC_NULL,"-da_view_draw",&flg1);
1751:   if (flg1) {DAView(da,PETSC_VIEWER_DRAW_(da->comm));}
1752:   PetscOptionsHasName(PETSC_NULL,"-help",&flg1);
1753:   if (flg1) {DAPrintHelp(da);}
1754:   PetscPublishAll(da);

1756: #if defined(PETSC_HAVE_AMS)
1757:   PetscObjectComposeFunctionDynamic((PetscObject)global,"AMSSetFieldBlock_C",
1758:          "AMSSetFieldBlock_DA",AMSSetFieldBlock_DA);
1759:   PetscObjectComposeFunctionDynamic((PetscObject)local,"AMSSetFieldBlock_C",
1760:          "AMSSetFieldBlock_DA",AMSSetFieldBlock_DA);
1761:   if (((PetscObject)global)->amem > -1) {
1762:     AMSSetFieldBlock_DA(((PetscObject)global)->amem,"values",global);
1763:   }
1764: #endif
1765:   VecSetOperation(global,VECOP_VIEW,(void(*)(void))VecView_MPI_DA);
1766:   VecSetOperation(global,VECOP_LOADINTOVECTOR,(void(*)(void))VecLoadIntoVector_Binary_DA);
1767:   return(0);
1768: }