Actual source code: da3.c

  1: /*$Id: da3.c,v 1.132 2001/03/28 19:42:42 balay 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"     /*I   "petscda.h"    I*/

 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:       double *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:     double     ymin = -1.0,ymax = (double)da->N;
 48:     double     xmin = -1.0,xmax = (double)((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 = (double)(da->N - 1);
 62:         for (xmin=(double)(k*(da->M+1)); xmin<(double)(da->M+(k*(da->M+1))); xmin++) {
 63:           PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);
 64:         }
 65: 
 66:         xmin = (double)(k*(da->M+1)); xmax = xmin + (double)(da->M - 1);
 67:         for (ymin=0; ymin<(double)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
181: .  -da_grid_y <ny> - number of grid points in y direction
182: .  -da_grid_z <nz> - number of grid points in z direction
183: -  -da_noao - do not compute natural to PETSc ordering object

185:    Level: beginner

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

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

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

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

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

223:   *inra = 0;

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

231:   PetscOptionsBegin(comm,PETSC_NULL,"3d DA Options","DA");
232:     PetscOptionsInt("-da_grid_x","Number of grid points in x direction","DACreate3d",M,&M,PETSC_NULL);
233:     PetscOptionsInt("-da_grid_y","Number of grid points in y direction","DACreate3d",N,&N,PETSC_NULL);
234:     PetscOptionsInt("-da_grid_z","Number of grid points in z direction","DACreate3d",P,&P,PETSC_NULL);
235:     PetscOptionsInt("-da_processors_x","Number of processors in x direction","DACreate3d",m,&m,PETSC_NULL);
236:     PetscOptionsInt("-da_processors_y","Number of processors in y direction","DACreate3d",n,&n,PETSC_NULL);
237:     PetscOptionsInt("-da_processors_z","Number of processors in z direction","DACreate3d",p,&p,PETSC_NULL);
238:   PetscOptionsEnd();

240:   PetscHeaderCreate(da,_p_DA,struct _DAOps,DA_COOKIE,0,"DA",comm,DADestroy,DAView);
241:   da->bops->publish           = DAPublish_Petsc;
242:   da->ops->createglobalvector = DACreateGlobalVector;
243:   da->ops->getinterpolation   = DAGetInterpolation;
244:   da->ops->getcoloring        = DAGetColoring;
245:   da->ops->refine             = DARefine;

247:   PetscLogObjectCreate(da);
248:   PetscLogObjectMemory(da,sizeof(struct _p_DA));
249:   da->dim        = 3;
250:   da->gtog1      = 0;
251:   PetscMalloc(dof*sizeof(char*),&da->fieldname);
252:   PetscMemzero(da->fieldname,dof*sizeof(char*));

254:   MPI_Comm_size(comm,&size);
255:   MPI_Comm_rank(comm,&rank);

257:   if (m != PETSC_DECIDE) {
258:     if (m < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %d",m);}
259:     else if (m > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %d %d",m,size);}
260:   }
261:   if (n != PETSC_DECIDE) {
262:     if (n < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %d",n);}
263:     else if (n > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %d %d",n,size);}
264:   }
265:   if (p != PETSC_DECIDE) {
266:     if (p < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %d",p);}
267:     else if (p > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %d %d",p,size);}
268:   }

270:   /* Partition the array among the processors */
271:   if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) {
272:     m = size/(n*p);
273:   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
274:     n = size/(m*p);
275:   } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
276:     p = size/(m*n);
277:   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
278:     /* try for squarish distribution */
279:     m = (int)(0.5 + sqrt(((double)M)*((double)size)/((double)N*p)));
280:     if (!m) m = 1;
281:     while (m > 0) {
282:       n = size/(m*p);
283:       if (m*n*p == size) break;
284:       m--;
285:     }
286:     if (!m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %d",p);
287:     if (M > N && m < n) {int _m = m; m = n; n = _m;}
288:   } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
289:     /* try for squarish distribution */
290:     m = (int)(0.5 + sqrt(((double)M)*((double)size)/((double)P*n)));
291:     if (!m) m = 1;
292:     while (m > 0) {
293:       p = size/(m*n);
294:       if (m*n*p == size) break;
295:       m--;
296:     }
297:     if (!m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %d",n);
298:     if (M > P && m < p) {int _m = m; m = p; p = _m;}
299:   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
300:     /* try for squarish distribution */
301:     n = (int)(0.5 + sqrt(((double)N)*((double)size)/((double)P*m)));
302:     if (!n) n = 1;
303:     while (n > 0) {
304:       p = size/(m*n);
305:       if (m*n*p == size) break;
306:       n--;
307:     }
308:     if (!n) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %d",n);
309:     if (N > P && n < p) {int _n = n; n = p; p = _n;}
310:   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
311:     /* try for squarish distribution */
312:     n = (int)(0.5 + pow(((double)N*N)*((double)size)/((double)P*M),1./3.));
313:     if (!n) n = 1;
314:     while (n > 0) {
315:       pm = size/n;
316:       if (n*pm == size) break;
317:       n--;
318:     }
319:     if (!n) n = 1;
320:     m = (int)(0.5 + sqrt(((double)M)*((double)size)/((double)P*n)));
321:     if (!m) m = 1;
322:     while (m > 0) {
323:       p = size/(m*n);
324:       if (m*n*p == size) break;
325:       m--;
326:     }
327:     if (M > P && m < p) {int _m = m; m = p; p = _m;}
328:   } else if (m*n*p != size) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");

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

335:   PetscOptionsHasName(PETSC_NULL,"-da_partition_nodes_at_end",&flg2);
336:   /* 
337:      Determine locally owned region 
338:      [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes 
339:   */
340:   if (lx) { /* user decided distribution */
341:     x  = lx[rank % m];
342:     xs = 0;
343:     for (i=0; i<(rank%m); i++) { xs += lx[i];}
344:     if (x < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column width is too thin for stencil! %d %d",x,s);
345:   } else if (flg2) {
346:     SETERRQ(PETSC_ERR_SUP,"-da_partition_nodes_at_end not supported");
347:   } else { /* Normal PETSc distribution */
348:     x = M/m + ((M % m) > (rank % m));
349:     if (x < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column width is too thin for stencil! %d %d",x,s);
350:     if ((M % m) > (rank % m)) { xs = (rank % m)*x; }
351:     else                      { xs = (M % m)*(x+1) + ((rank % m)-(M % m))*x; }
352:     PetscMalloc(m*sizeof(int),&lx);
353:     flx = lx;
354:     for (i=0; i<m; i++) {
355:       lx[i] = M/m + ((M % m) > (i % m));
356:     }
357:   }
358:   if (ly) { /* user decided distribution */
359:     y  = ly[(rank % (m*n))/m];
360:     if (y < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row width is too thin for stencil! %d %d",y,s);
361:     ys = 0;
362:     for (i=0; i<(rank % (m*n))/m; i++) { ys += ly[i];}
363:   } else if (flg2) {
364:     SETERRQ(PETSC_ERR_SUP,"-da_partition_nodes_at_end not supported");
365:   } else { /* Normal PETSc distribution */
366:     y = N/n + ((N % n) > ((rank % (m*n)) /m));
367:     if (y < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row width is too thin for stencil! %d %d",y,s);
368:     if ((N % n) > ((rank % (m*n)) /m)) {ys = ((rank % (m*n))/m)*y;}
369:     else                               {ys = (N % n)*(y+1) + (((rank % (m*n))/m)-(N % n))*y;}
370:     PetscMalloc(n*sizeof(int),&ly);
371:     fly = ly;
372:     for (i=0; i<n; i++) {
373:       ly[i] = N/n + ((N % n) > (i % n));
374:     }
375:   }
376:   if (lz) { /* user decided distribution */
377:     z  = lz[rank/(m*n)];
378:     if (z < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Plane width is too thin for stencil! %d %d",z,s);
379:     zs = 0;
380:     for (i=0; i<(rank/(m*n)); i++) { zs += lz[i];}
381:   } else if (flg2) {
382:     SETERRQ(PETSC_ERR_SUP,"-da_partition_nodes_at_end not supported");
383:   } else { /* Normal PETSc distribution */
384:     z = P/p + ((P % p) > (rank / (m*n)));
385:     if (z < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Plane width is too thin for stencil! %d %d",z,s);
386:     if ((P % p) > (rank / (m*n))) {zs = (rank/(m*n))*z;}
387:     else                          {zs = (P % p)*(z+1) + ((rank/(m*n))-(P % p))*z;}
388:     PetscMalloc(p*sizeof(int),&lz);
389:     flz = lz;
390:     for (i=0; i<p; i++) {
391:       lz[i] = P/p + ((P % p) > (i % p));
392:     }
393:   }
394:   ye = ys + y;
395:   xe = xs + x;
396:   ze = zs + z;

398:   /* determine ghost region */
399:   /* Assume No Periodicity */
400:   if (xs-s > 0) Xs = xs - s; else Xs = 0;
401:   if (ys-s > 0) Ys = ys - s; else Ys = 0;
402:   if (zs-s > 0) Zs = zs - s; else Zs = 0;
403:   if (xe+s <= M) Xe = xe + s; else Xe = M;
404:   if (ye+s <= N) Ye = ye + s; else Ye = N;
405:   if (ze+s <= P) Ze = ze + s; else Ze = P;

407:   /* X Periodic */
408:   if (DAXPeriodic(wrap)){
409:     Xs = xs - s;
410:     Xe = xe + s;
411:   }

413:   /* Y Periodic */
414:   if (DAYPeriodic(wrap)){
415:     Ys = ys - s;
416:     Ye = ye + s;
417:   }

419:   /* Z Periodic */
420:   if (DAZPeriodic(wrap)){
421:     Zs = zs - s;
422:     Ze = ze + s;
423:   }

425:   /* Resize all X parameters to reflect w */
426:   x   *= dof;
427:   xs  *= dof;
428:   xe  *= dof;
429:   Xs  *= dof;
430:   Xe  *= dof;
431:   s_x  = s*dof;
432:   s_y  = s;
433:   s_z  = s;

435:   /* determine starting point of each processor */
436:   nn = x*y*z;
437:   PetscMalloc((2*size+1)*sizeof(int),&bases);
438:   ldims = (int*)(bases+size+1);
439:   MPI_Allgather(&nn,1,MPI_INT,ldims,1,MPI_INT,comm);
440:   bases[0] = 0;
441:   for (i=1; i<=size; i++) {
442:     bases[i] = ldims[i-1];
443:   }
444:   for (i=1; i<=size; i++) {
445:     bases[i] += bases[i-1];
446:   }

448:   /* allocate the base parallel and sequential vectors */
449:   VecCreateMPI(comm,x*y*z,PETSC_DECIDE,&global);
450:   VecSetBlockSize(global,dof);
451:   VecCreateSeq(MPI_COMM_SELF,(Xe-Xs)*(Ye-Ys)*(Ze-Zs),&local);
452:   VecSetBlockSize(local,dof);

454:   /* generate appropriate vector scatters */
455:   /* local to global inserts non-ghost point region into global */
456:   VecGetOwnershipRange(global,&start,&end);
457:   ISCreateStride(comm,x*y*z,start,1,&to);

459:   left   = xs - Xs;
460:   bottom = ys - Ys; top = bottom + y;
461:   down   = zs - Zs; up  = down + z;
462:   count  = x*(top-bottom)*(up-down);
463:   PetscMalloc(count*sizeof(int),&idx);
464:   count  = 0;
465:   for (i=down; i<up; i++) {
466:     for (j=bottom; j<top; j++) {
467:       for (k=0; k<x; k++) {
468:         idx[count++] = (left+j*(Xe-Xs))+i*(Xe-Xs)*(Ye-Ys) + k;
469:       }
470:     }
471:   }
472:   ISCreateGeneral(comm,count,idx,&from);
473:   PetscFree(idx);

475:   VecScatterCreate(local,from,global,to,&ltog);
476:   PetscLogObjectParent(da,to);
477:   PetscLogObjectParent(da,from);
478:   PetscLogObjectParent(da,ltog);
479:   ISDestroy(from);
480:   ISDestroy(to);

482:   /* global to local must include ghost points */
483:   if (stencil_type == DA_STENCIL_BOX) {
484:     ISCreateStride(comm,(Xe-Xs)*(Ye-Ys)*(Ze-Zs),0,1,&to);
485:   } else {
486:     /* This is way ugly! We need to list the funny cross type region */
487:     /* the bottom chunck */
488:     left   = xs - Xs;
489:     bottom = ys - Ys; top = bottom + y;
490:     down   = zs - Zs;   up  = down + z;
491:     count  = down*(top-bottom)*x +
492:              (up-down)*(bottom*x  + (top-bottom)*(Xe-Xs) + (Ye-Ys-top)*x) +
493:              (Ze-Zs-up)*(top-bottom)*x;
494:     PetscMalloc(count*sizeof(int),&idx);
495:     count  = 0;
496:     for (i=0; i<down; i++) {
497:       for (j=bottom; j<top; j++) {
498:         for (k=0; k<x; k++) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
499:       }
500:     }
501:     /* the middle piece */
502:     for (i=down; i<up; i++) {
503:       /* front */
504:       for (j=0; j<bottom; j++) {
505:         for (k=0; k<x; k++) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
506:       }
507:       /* middle */
508:       for (j=bottom; j<top; j++) {
509:         for (k=0; k<Xe-Xs; k++) idx[count++] = j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
510:       }
511:       /* back */
512:       for (j=top; j<Ye-Ys; j++) {
513:         for (k=0; k<x; k++) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
514:       }
515:     }
516:     /* the top piece */
517:     for (i=up; i<Ze-Zs; i++) {
518:       for (j=bottom; j<top; j++) {
519:         for (k=0; k<x; k++) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
520:       }
521:     }
522:     ISCreateGeneral(comm,count,idx,&to);
523:     PetscFree(idx);
524:   }

526:   /* determine who lies on each side of use stored in    n24 n25 n26
527:                                                          n21 n22 n23
528:                                                          n18 n19 n20

530:                                                          n15 n16 n17
531:                                                          n12     n14
532:                                                          n9  n10 n11

534:                                                          n6  n7  n8
535:                                                          n3  n4  n5
536:                                                          n0  n1  n2
537:   */
538: 
539:   /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
540: 
541:   /* Assume Nodes are Internal to the Cube */
542: 
543:   n0  = rank - m*n - m - 1;
544:   n1  = rank - m*n - m;
545:   n2  = rank - m*n - m + 1;
546:   n3  = rank - m*n -1;
547:   n4  = rank - m*n;
548:   n5  = rank - m*n + 1;
549:   n6  = rank - m*n + m - 1;
550:   n7  = rank - m*n + m;
551:   n8  = rank - m*n + m + 1;

553:   n9  = rank - m - 1;
554:   n10 = rank - m;
555:   n11 = rank - m + 1;
556:   n12 = rank - 1;
557:   n14 = rank + 1;
558:   n15 = rank + m - 1;
559:   n16 = rank + m;
560:   n17 = rank + m + 1;

562:   n18 = rank + m*n - m - 1;
563:   n19 = rank + m*n - m;
564:   n20 = rank + m*n - m + 1;
565:   n21 = rank + m*n - 1;
566:   n22 = rank + m*n;
567:   n23 = rank + m*n + 1;
568:   n24 = rank + m*n + m - 1;
569:   n25 = rank + m*n + m;
570:   n26 = rank + m*n + m + 1;

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

574:   if (xs == 0) { /* First assume not corner or edge */
575:     n0  = rank       -1 - (m*n);
576:     n3  = rank + m   -1 - (m*n);
577:     n6  = rank + 2*m -1 - (m*n);
578:     n9  = rank       -1;
579:     n12 = rank + m   -1;
580:     n15 = rank + 2*m -1;
581:     n18 = rank       -1 + (m*n);
582:     n21 = rank + m   -1 + (m*n);
583:     n24 = rank + 2*m -1 + (m*n);
584:    }

586:   if (xe == M*dof) { /* First assume not corner or edge */
587:     n2  = rank -2*m +1 - (m*n);
588:     n5  = rank - m  +1 - (m*n);
589:     n8  = rank      +1 - (m*n);
590:     n11 = rank -2*m +1;
591:     n14 = rank - m  +1;
592:     n17 = rank      +1;
593:     n20 = rank -2*m +1 + (m*n);
594:     n23 = rank - m  +1 + (m*n);
595:     n26 = rank      +1 + (m*n);
596:   }

598:   if (ys==0) { /* First assume not corner or edge */
599:     n0  = rank + m * (n-1) -1 - (m*n);
600:     n1  = rank + m * (n-1)    - (m*n);
601:     n2  = rank + m * (n-1) +1 - (m*n);
602:     n9  = rank + m * (n-1) -1;
603:     n10 = rank + m * (n-1);
604:     n11 = rank + m * (n-1) +1;
605:     n18 = rank + m * (n-1) -1 + (m*n);
606:     n19 = rank + m * (n-1)    + (m*n);
607:     n20 = rank + m * (n-1) +1 + (m*n);
608:   }

610:   if (ye == N) { /* First assume not corner or edge */
611:     n6  = rank - m * (n-1) -1 - (m*n);
612:     n7  = rank - m * (n-1)    - (m*n);
613:     n8  = rank - m * (n-1) +1 - (m*n);
614:     n15 = rank - m * (n-1) -1;
615:     n16 = rank - m * (n-1);
616:     n17 = rank - m * (n-1) +1;
617:     n24 = rank - m * (n-1) -1 + (m*n);
618:     n25 = rank - m * (n-1)    + (m*n);
619:     n26 = rank - m * (n-1) +1 + (m*n);
620:   }
621: 
622:   if (zs == 0) { /* First assume not corner or edge */
623:     n0 = size - (m*n) + rank - m - 1;
624:     n1 = size - (m*n) + rank - m;
625:     n2 = size - (m*n) + rank - m + 1;
626:     n3 = size - (m*n) + rank - 1;
627:     n4 = size - (m*n) + rank;
628:     n5 = size - (m*n) + rank + 1;
629:     n6 = size - (m*n) + rank + m - 1;
630:     n7 = size - (m*n) + rank + m ;
631:     n8 = size - (m*n) + rank + m + 1;
632:   }

634:   if (ze == P) { /* First assume not corner or edge */
635:     n18 = (m*n) - (size-rank) - m - 1;
636:     n19 = (m*n) - (size-rank) - m;
637:     n20 = (m*n) - (size-rank) - m + 1;
638:     n21 = (m*n) - (size-rank) - 1;
639:     n22 = (m*n) - (size-rank);
640:     n23 = (m*n) - (size-rank) + 1;
641:     n24 = (m*n) - (size-rank) + m - 1;
642:     n25 = (m*n) - (size-rank) + m;
643:     n26 = (m*n) - (size-rank) + m + 1;
644:   }

646:   if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */
647:     n0 = size - m*n + rank + m-1 - m;
648:     n3 = size - m*n + rank + m-1;
649:     n6 = size - m*n + rank + m-1 + m;
650:   }
651: 
652:   if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */
653:     n18 = m*n - (size - rank) + m-1 - m;
654:     n21 = m*n - (size - rank) + m-1;
655:     n24 = m*n - (size - rank) + m-1 + m;
656:   }

658:   if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */
659:     n0  = rank + m*n -1 - m*n;
660:     n9  = rank + m*n -1;
661:     n18 = rank + m*n -1 + m*n;
662:   }

664:   if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */
665:     n6  = rank - m*(n-1) + m-1 - m*n;
666:     n15 = rank - m*(n-1) + m-1;
667:     n24 = rank - m*(n-1) + m-1 + m*n;
668:   }

670:   if ((xe==M*dof) && (zs==0)) { /* Assume an edge, not corner */
671:     n2 = size - (m*n-rank) - (m-1) - m;
672:     n5 = size - (m*n-rank) - (m-1);
673:     n8 = size - (m*n-rank) - (m-1) + m;
674:   }

676:   if ((xe==M*dof) && (ze==P)) { /* Assume an edge, not corner */
677:     n20 = m*n - (size - rank) - (m-1) - m;
678:     n23 = m*n - (size - rank) - (m-1);
679:     n26 = m*n - (size - rank) - (m-1) + m;
680:   }

682:   if ((xe==M*dof) && (ys==0)) { /* Assume an edge, not corner */
683:     n2  = rank + m*(n-1) - (m-1) - m*n;
684:     n11 = rank + m*(n-1) - (m-1);
685:     n20 = rank + m*(n-1) - (m-1) + m*n;
686:   }

688:   if ((xe==M*dof) && (ye==N)) { /* Assume an edge, not corner */
689:     n8  = rank - m*n +1 - m*n;
690:     n17 = rank - m*n +1;
691:     n26 = rank - m*n +1 + m*n;
692:   }

694:   if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */
695:     n0 = size - m + rank -1;
696:     n1 = size - m + rank;
697:     n2 = size - m + rank +1;
698:   }

700:   if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */
701:     n18 = m*n - (size - rank) + m*(n-1) -1;
702:     n19 = m*n - (size - rank) + m*(n-1);
703:     n20 = m*n - (size - rank) + m*(n-1) +1;
704:   }

706:   if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */
707:     n6 = size - (m*n-rank) - m * (n-1) -1;
708:     n7 = size - (m*n-rank) - m * (n-1);
709:     n8 = size - (m*n-rank) - m * (n-1) +1;
710:   }

712:   if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */
713:     n24 = rank - (size-m) -1;
714:     n25 = rank - (size-m);
715:     n26 = rank - (size-m) +1;
716:   }

718:   /* Check for Corners */
719:   if ((xs==0)   && (ys==0) && (zs==0)) { n0  = size -1;}
720:   if ((xs==0)   && (ys==0) && (ze==P)) { n18 = m*n-1;}
721:   if ((xs==0)   && (ye==N) && (zs==0)) { n6  = (size-1)-m*(n-1);}
722:   if ((xs==0)   && (ye==N) && (ze==P)) { n24 = m-1;}
723:   if ((xe==M*dof) && (ys==0) && (zs==0)) { n2  = size-m;}
724:   if ((xe==M*dof) && (ys==0) && (ze==P)) { n20 = m*n-m;}
725:   if ((xe==M*dof) && (ye==N) && (zs==0)) { n8  = size-m*n;}
726:   if ((xe==M*dof) && (ye==N) && (ze==P)) { n26 = 0;}

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

730:   /* If not X periodic */
731:   if ((wrap != DA_XPERIODIC)  && (wrap != DA_XYPERIODIC) &&
732:      (wrap != DA_XZPERIODIC) && (wrap != DA_XYZPERIODIC)) {
733:     if (xs==0)   {n0  = n3  = n6  = n9  = n12 = n15 = n18 = n21 = n24 = -2;}
734:     if (xe==M*dof) {n2  = n5  = n8  = n11 = n14 = n17 = n20 = n23 = n26 = -2;}
735:   }

737:   /* If not Y periodic */
738:   if ((wrap != DA_YPERIODIC)  && (wrap != DA_XYPERIODIC) &&
739:       (wrap != DA_YZPERIODIC) && (wrap != DA_XYZPERIODIC)) {
740:     if (ys==0)   {n0  = n1  = n2  = n9  = n10 = n11 = n18 = n19 = n20 = -2;}
741:     if (ye==N)   {n6  = n7  = n8  = n15 = n16 = n17 = n24 = n25 = n26 = -2;}
742:   }

744:   /* If not Z periodic */
745:   if ((wrap != DA_ZPERIODIC)  && (wrap != DA_XZPERIODIC) &&
746:       (wrap != DA_YZPERIODIC) && (wrap != DA_XYZPERIODIC)) {
747:     if (zs==0)   {n0  = n1  = n2  = n3  = n4  = n5  = n6  = n7  = n8  = -2;}
748:     if (ze==P)   {n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;}
749:   }

751:   /* If star stencil then delete the corner neighbors */
752:   if (stencil_type == DA_STENCIL_STAR) {
753:      /* save information about corner neighbors */
754:      sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7;
755:      sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18;
756:      sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25;
757:      sn26 = n26;
758:      n0  = n1  = n2  = n3  = n5  = n6  = n7  = n8  = n9  = n11 =
759:      n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1;
760:   }


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

766:   nn = 0;

768:   /* Bottom Level */
769:   for (k=0; k<s_z; k++) {
770:     for (i=1; i<=s_y; i++) {
771:       if (n0 >= 0) { /* left below */
772:         x_t = lx[n0 % m]*dof;
773:         y_t = ly[(n0 % (m*n))/m];
774:         z_t = lz[n0 / (m*n)];
775:         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;
776:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
777:       }
778:       if (n1 >= 0) { /* directly below */
779:         x_t = x;
780:         y_t = ly[(n1 % (m*n))/m];
781:         z_t = lz[n1 / (m*n)];
782:         s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
783:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
784:       }
785:       if (n2 >= 0) { /* right below */
786:         x_t = lx[n2 % m]*dof;
787:         y_t = ly[(n2 % (m*n))/m];
788:         z_t = lz[n2 / (m*n)];
789:         s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
790:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
791:       }
792:     }

794:     for (i=0; i<y; i++) {
795:       if (n3 >= 0) { /* directly left */
796:         x_t = lx[n3 % m]*dof;
797:         y_t = y;
798:         z_t = lz[n3 / (m*n)];
799:         s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
800:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
801:       }

803:       if (n4 >= 0) { /* middle */
804:         x_t = x;
805:         y_t = y;
806:         z_t = lz[n4 / (m*n)];
807:         s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
808:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
809:       }

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

820:     for (i=1; i<=s_y; i++) {
821:       if (n6 >= 0) { /* left above */
822:         x_t = lx[n6 % m]*dof;
823:         y_t = ly[(n6 % (m*n))/m];
824:         z_t = lz[n6 / (m*n)];
825:         s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
826:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
827:       }
828:       if (n7 >= 0) { /* directly above */
829:         x_t = x;
830:         y_t = ly[(n7 % (m*n))/m];
831:         z_t = lz[n7 / (m*n)];
832:         s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
833:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
834:       }
835:       if (n8 >= 0) { /* right above */
836:         x_t = lx[n8 % m]*dof;
837:         y_t = ly[(n8 % (m*n))/m];
838:         z_t = lz[n8 / (m*n)];
839:         s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
840:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
841:       }
842:     }
843:   }

845:   /* Middle Level */
846:   for (k=0; k<z; k++) {
847:     for (i=1; i<=s_y; i++) {
848:       if (n9 >= 0) { /* left below */
849:         x_t = lx[n9 % m]*dof;
850:         y_t = ly[(n9 % (m*n))/m];
851:         /* z_t = z; */
852:         s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
853:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
854:       }
855:       if (n10 >= 0) { /* directly below */
856:         x_t = x;
857:         y_t = ly[(n10 % (m*n))/m];
858:         /* z_t = z; */
859:         s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
860:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
861:       }
862:       if (n11 >= 0) { /* right below */
863:         x_t = lx[n11 % m]*dof;
864:         y_t = ly[(n11 % (m*n))/m];
865:         /* z_t = z; */
866:         s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
867:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
868:       }
869:     }

871:     for (i=0; i<y; i++) {
872:       if (n12 >= 0) { /* directly left */
873:         x_t = lx[n12 % m]*dof;
874:         y_t = y;
875:         /* z_t = z; */
876:         s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
877:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
878:       }

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

884:       if (n14 >= 0) { /* directly right */
885:         x_t = lx[n14 % m]*dof;
886:         y_t = y;
887:         /* z_t = z; */
888:         s_t = bases[n14] + i*x_t + k*x_t*y_t;
889:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
890:       }
891:     }

893:     for (i=1; i<=s_y; i++) {
894:       if (n15 >= 0) { /* left above */
895:         x_t = lx[n15 % m]*dof;
896:         y_t = ly[(n15 % (m*n))/m];
897:         /* z_t = z; */
898:         s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
899:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
900:       }
901:       if (n16 >= 0) { /* directly above */
902:         x_t = x;
903:         y_t = ly[(n16 % (m*n))/m];
904:         /* z_t = z; */
905:         s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
906:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
907:       }
908:       if (n17 >= 0) { /* right above */
909:         x_t = lx[n17 % m]*dof;
910:         y_t = ly[(n17 % (m*n))/m];
911:         /* z_t = z; */
912:         s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
913:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
914:       }
915:     }
916:   }
917: 
918:   /* Upper Level */
919:   for (k=0; k<s_z; k++) {
920:     for (i=1; i<=s_y; i++) {
921:       if (n18 >= 0) { /* left below */
922:         x_t = lx[n18 % m]*dof;
923:         y_t = ly[(n18 % (m*n))/m];
924:         /* z_t = lz[n18 / (m*n)]; */
925:         s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
926:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
927:       }
928:       if (n19 >= 0) { /* directly below */
929:         x_t = x;
930:         y_t = ly[(n19 % (m*n))/m];
931:         /* z_t = lz[n19 / (m*n)]; */
932:         s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
933:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
934:       }
935:       if (n20 >= 0) { /* right below */
936:         x_t = lx[n20 % m]*dof;
937:         y_t = ly[(n20 % (m*n))/m];
938:         /* z_t = lz[n20 / (m*n)]; */
939:         s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
940:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
941:       }
942:     }

944:     for (i=0; i<y; i++) {
945:       if (n21 >= 0) { /* directly left */
946:         x_t = lx[n21 % m]*dof;
947:         y_t = y;
948:         /* z_t = lz[n21 / (m*n)]; */
949:         s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
950:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
951:       }

953:       if (n22 >= 0) { /* middle */
954:         x_t = x;
955:         y_t = y;
956:         /* z_t = lz[n22 / (m*n)]; */
957:         s_t = bases[n22] + i*x_t + k*x_t*y_t;
958:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
959:       }

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

970:     for (i=1; i<=s_y; i++) {
971:       if (n24 >= 0) { /* left above */
972:         x_t = lx[n24 % m]*dof;
973:         y_t = ly[(n24 % (m*n))/m];
974:         /* z_t = lz[n24 / (m*n)]; */
975:         s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
976:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
977:       }
978:       if (n25 >= 0) { /* directly above */
979:         x_t = x;
980:         y_t = ly[(n25 % (m*n))/m];
981:         /* z_t = lz[n25 / (m*n)]; */
982:         s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
983:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
984:       }
985:       if (n26 >= 0) { /* right above */
986:         x_t = lx[n26 % m]*dof;
987:         y_t = ly[(n26 % (m*n))/m];
988:         /* z_t = lz[n26 / (m*n)]; */
989:         s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
990:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
991:       }
992:     }
993:   }
994:   base = bases[rank];
995:   ISCreateGeneral(comm,nn,idx,&from);
996:   VecScatterCreate(global,from,local,to,&gtol);
997:   PetscLogObjectParent(da,gtol);
998:   PetscLogObjectParent(da,to);
999:   PetscLogObjectParent(da,from);
1000:   ISDestroy(to);
1001:   ISDestroy(from);
1002:   da->stencil_type = stencil_type;
1003:   da->M  = M;  da->N  = N; da->P = P;
1004:   da->m  = m;  da->n  = n; da->p = p;
1005:   da->w  = dof;  da->s  = s;
1006:   da->xs = xs; da->xe = xe; da->ys = ys; da->ye = ye; da->zs = zs; da->ze = ze;
1007:   da->Xs = Xs; da->Xe = Xe; da->Ys = Ys; da->Ye = Ye; da->Zs = Zs; da->Ze = Ze;

1009:   PetscLogObjectParent(da,global);
1010:   PetscLogObjectParent(da,local);

1012:   if (stencil_type == DA_STENCIL_STAR) {
1013:     /*
1014:         Recompute the local to global mappings, this time keeping the 
1015:       information about the cross corner processor numbers.
1016:     */
1017:     n0  = sn0;  n1  = sn1;  n2  = sn2;  n3  = sn3;  n5  = sn5;  n6  = sn6; n7 = sn7;
1018:     n8  = sn8;  n9  = sn9;  n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18;
1019:     n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25;
1020:     n26 = sn26;

1022:     nn = 0;

1024:     /* Bottom Level */
1025:     for (k=0; k<s_z; k++) {
1026:       for (i=1; i<=s_y; i++) {
1027:         if (n0 >= 0) { /* left below */
1028:           x_t = lx[n0 % m]*dof;
1029:           y_t = ly[(n0 % (m*n))/m];
1030:           z_t = lz[n0 / (m*n)];
1031:           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;
1032:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1033:         }
1034:         if (n1 >= 0) { /* directly below */
1035:           x_t = x;
1036:           y_t = ly[(n1 % (m*n))/m];
1037:           z_t = lz[n1 / (m*n)];
1038:           s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
1039:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1040:         }
1041:         if (n2 >= 0) { /* right below */
1042:           x_t = lx[n2 % m]*dof;
1043:           y_t = ly[(n2 % (m*n))/m];
1044:           z_t = lz[n2 / (m*n)];
1045:           s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
1046:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1047:         }
1048:       }

1050:       for (i=0; i<y; i++) {
1051:         if (n3 >= 0) { /* directly left */
1052:           x_t = lx[n3 % m]*dof;
1053:           y_t = y;
1054:           z_t = lz[n3 / (m*n)];
1055:           s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1056:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1057:         }

1059:         if (n4 >= 0) { /* middle */
1060:           x_t = x;
1061:           y_t = y;
1062:           z_t = lz[n4 / (m*n)];
1063:           s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1064:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1065:         }

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

1076:       for (i=1; i<=s_y; i++) {
1077:         if (n6 >= 0) { /* left above */
1078:           x_t = lx[n6 % m]*dof;
1079:           y_t = ly[(n6 % (m*n))/m];
1080:           z_t = lz[n6 / (m*n)];
1081:           s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1082:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1083:         }
1084:         if (n7 >= 0) { /* directly above */
1085:           x_t = x;
1086:           y_t = ly[(n7 % (m*n))/m];
1087:           z_t = lz[n7 / (m*n)];
1088:           s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1089:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1090:         }
1091:         if (n8 >= 0) { /* right above */
1092:           x_t = lx[n8 % m]*dof;
1093:           y_t = ly[(n8 % (m*n))/m];
1094:           z_t = lz[n8 / (m*n)];
1095:           s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1096:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1097:         }
1098:       }
1099:     }

1101:     /* Middle Level */
1102:     for (k=0; k<z; k++) {
1103:       for (i=1; i<=s_y; i++) {
1104:         if (n9 >= 0) { /* left below */
1105:           x_t = lx[n9 % m]*dof;
1106:           y_t = ly[(n9 % (m*n))/m];
1107:           /* z_t = z; */
1108:           s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
1109:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1110:         }
1111:         if (n10 >= 0) { /* directly below */
1112:           x_t = x;
1113:           y_t = ly[(n10 % (m*n))/m];
1114:           /* z_t = z; */
1115:           s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1116:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1117:         }
1118:         if (n11 >= 0) { /* right below */
1119:           x_t = lx[n11 % m]*dof;
1120:           y_t = ly[(n11 % (m*n))/m];
1121:           /* z_t = z; */
1122:           s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1123:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1124:         }
1125:       }

1127:       for (i=0; i<y; i++) {
1128:         if (n12 >= 0) { /* directly left */
1129:           x_t = lx[n12 % m]*dof;
1130:           y_t = y;
1131:           /* z_t = z; */
1132:           s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
1133:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1134:         }

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

1140:         if (n14 >= 0) { /* directly right */
1141:           x_t = lx[n14 % m]*dof;
1142:           y_t = y;
1143:           /* z_t = z; */
1144:           s_t = bases[n14] + i*x_t + k*x_t*y_t;
1145:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1146:         }
1147:       }

1149:       for (i=1; i<=s_y; i++) {
1150:         if (n15 >= 0) { /* left above */
1151:           x_t = lx[n15 % m]*dof;
1152:           y_t = ly[(n15 % (m*n))/m];
1153:           /* z_t = z; */
1154:           s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
1155:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1156:         }
1157:         if (n16 >= 0) { /* directly above */
1158:           x_t = x;
1159:           y_t = ly[(n16 % (m*n))/m];
1160:           /* z_t = z; */
1161:           s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
1162:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1163:         }
1164:         if (n17 >= 0) { /* right above */
1165:           x_t = lx[n17 % m]*dof;
1166:           y_t = ly[(n17 % (m*n))/m];
1167:           /* z_t = z; */
1168:           s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
1169:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1170:         }
1171:       }
1172:     }
1173: 
1174:     /* Upper Level */
1175:     for (k=0; k<s_z; k++) {
1176:       for (i=1; i<=s_y; i++) {
1177:         if (n18 >= 0) { /* left below */
1178:           x_t = lx[n18 % m]*dof;
1179:           y_t = ly[(n18 % (m*n))/m];
1180:           /* z_t = lz[n18 / (m*n)]; */
1181:           s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
1182:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1183:         }
1184:         if (n19 >= 0) { /* directly below */
1185:           x_t = x;
1186:           y_t = ly[(n19 % (m*n))/m];
1187:           /* z_t = lz[n19 / (m*n)]; */
1188:           s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1189:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1190:         }
1191:         if (n20 >= 0) { /* right below */
1192:           x_t = lx[n20 % m]*dof;
1193:           y_t = ly[(n20 % (m*n))/m];
1194:           /* z_t = lz[n20 / (m*n)]; */
1195:           s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1196:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1197:         }
1198:       }

1200:       for (i=0; i<y; i++) {
1201:         if (n21 >= 0) { /* directly left */
1202:           x_t = lx[n21 % m]*dof;
1203:           y_t = y;
1204:           /* z_t = lz[n21 / (m*n)]; */
1205:           s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
1206:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1207:         }

1209:         if (n22 >= 0) { /* middle */
1210:           x_t = x;
1211:           y_t = y;
1212:           /* z_t = lz[n22 / (m*n)]; */
1213:           s_t = bases[n22] + i*x_t + k*x_t*y_t;
1214:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1215:         }

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

1226:       for (i=1; i<=s_y; i++) {
1227:         if (n24 >= 0) { /* left above */
1228:           x_t = lx[n24 % m]*dof;
1229:           y_t = ly[(n24 % (m*n))/m];
1230:           /* z_t = lz[n24 / (m*n)]; */
1231:           s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
1232:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1233:         }
1234:         if (n25 >= 0) { /* directly above */
1235:           x_t = x;
1236:           y_t = ly[(n25 % (m*n))/m];
1237:           /* z_t = lz[n25 / (m*n)]; */
1238:           s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
1239:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1240:         }
1241:         if (n26 >= 0) { /* right above */
1242:           x_t = lx[n26 % m]*dof;
1243:           y_t = ly[(n26 % (m*n))/m];
1244:           /* z_t = lz[n26 / (m*n)]; */
1245:           s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
1246:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1247:         }
1248:       }
1249:     }
1250:   }
1251:   da->global    = global;
1252:   da->local     = local;
1253:   da->gtol      = gtol;
1254:   da->ltog      = ltog;
1255:   da->idx       = idx;
1256:   da->Nl        = nn;
1257:   da->base      = base;
1258:   da->ops->view = DAView_3d;
1259:   da->wrap      = wrap;
1260:   *inra = da;

1262:   /* 
1263:      Set the local to global ordering in the global vector, this allows use
1264:      of VecSetValuesLocal().
1265:   */
1266:   ierr  = ISLocalToGlobalMappingCreate(comm,nn,idx,&da->ltogmap);
1267:   ierr  = VecSetLocalToGlobalMapping(da->global,da->ltogmap);
1268:   ISLocalToGlobalMappingBlock(da->ltogmap,da->w,&da->ltogmapb);
1269:   VecSetLocalToGlobalMappingBlock(da->global,da->ltogmapb);
1270:   PetscLogObjectParent(da,da->ltogmap);

1272:   /* redo idx to include "missing" ghost points */
1273:   /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
1274: 
1275:   /* Assume Nodes are Internal to the Cube */
1276: 
1277:   n0  = rank - m*n - m - 1;
1278:   n1  = rank - m*n - m;
1279:   n2  = rank - m*n - m + 1;
1280:   n3  = rank - m*n -1;
1281:   n4  = rank - m*n;
1282:   n5  = rank - m*n + 1;
1283:   n6  = rank - m*n + m - 1;
1284:   n7  = rank - m*n + m;
1285:   n8  = rank - m*n + m + 1;

1287:   n9  = rank - m - 1;
1288:   n10 = rank - m;
1289:   n11 = rank - m + 1;
1290:   n12 = rank - 1;
1291:   n14 = rank + 1;
1292:   n15 = rank + m - 1;
1293:   n16 = rank + m;
1294:   n17 = rank + m + 1;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1482:   nn = 0;

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

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

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

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

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

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

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

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

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

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

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

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

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

1686:     for (i=1; i<=s_y; i++) {
1687:       if (n24 >= 0) { /* left above */
1688:         x_t = lx[n24 % m]*dof;
1689:         y_t = ly[(n24 % (m*n))/m];
1690:         /* z_t = lz[n24 / (m*n)]; */
1691:         s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
1692:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1693:       }
1694:       if (n25 >= 0) { /* directly above */
1695:         x_t = x;
1696:         y_t = ly[(n25 % (m*n))/m];
1697:         /* z_t = lz[n25 / (m*n)]; */
1698:         s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
1699:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1700:       }
1701:       if (n26 >= 0) { /* right above */
1702:         x_t = lx[n26 % m]*dof;
1703:         y_t = ly[(n26 % (m*n))/m];
1704:         /* z_t = lz[n26 / (m*n)]; */
1705:         s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
1706:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1707:       }
1708:     }
1709:   }

1711:   /* construct the local to local scatter context */
1712:   /* 
1713:       We simply remap the values in the from part of 
1714:     global to local to read from an array with the ghost values 
1715:     rather then from the plan array.
1716:   */
1717:   VecScatterCopy(gtol,&da->ltol);
1718:   PetscLogObjectParent(da,da->ltol);
1719:   left   = xs - Xs;
1720:   bottom = ys - Ys; top = bottom + y;
1721:   down   = zs - Zs; up  = down + z;
1722:   count  = x*(top-bottom)*(up-down);
1723:   PetscMalloc(count*sizeof(int),&idx);
1724:   count  = 0;
1725:   for (i=down; i<up; i++) {
1726:     for (j=bottom; j<top; j++) {
1727:       for (k=0; k<x; k++) {
1728:         idx[count++] = (left+j*(Xe-Xs))+i*(Xe-Xs)*(Ye-Ys) + k;
1729:       }
1730:     }
1731:   }
1732:   VecScatterRemap(da->ltol,idx,PETSC_NULL);
1733:   PetscFree(idx);


1736:   /* 
1737:      Build the natural ordering to PETSc ordering mappings.
1738:   */
1739:   PetscOptionsHasName(PETSC_NULL,"-da_noao",&flg1);
1740:   if (!flg1) {
1741:     IS  ispetsc,isnatural;
1742:     int *lidx,lict = 0;
1743:     int Nlocal = (da->xe-da->xs)*(da->ye-da->ys)*(da->ze-da->zs);

1745:     ISCreateStride(comm,Nlocal,da->base,1,&ispetsc);

1747:     PetscMalloc(Nlocal*sizeof(int),&lidx);
1748:     for (k=zs; k<ze; k++) {
1749:       for (j=ys; j<ye; j++) {
1750:         for (i=xs; i<xe; i++) {
1751:           lidx[lict++] = i + j*M*dof + k*M*N*dof;
1752:         }
1753:       }
1754:     }
1755:     ISCreateGeneral(comm,Nlocal,lidx,&isnatural);
1756:     PetscFree(lidx);

1758:     AOCreateBasicIS(isnatural,ispetsc,&da->ao);
1759:     PetscLogObjectParent(da,da->ao);
1760:     ISDestroy(ispetsc);
1761:     ISDestroy(isnatural);
1762:   } else {
1763:     da->ao = PETSC_NULL;
1764:   }

1766:   if (!flx) {
1767:     PetscMalloc(m*sizeof(int),&flx);
1768:     PetscMemcpy(flx,lx,m*sizeof(int));
1769:   }
1770:   if (!fly) {
1771:     PetscMalloc(n*sizeof(int),&fly);
1772:     PetscMemcpy(fly,ly,n*sizeof(int));
1773:   }
1774:   if (!flz) {
1775:     PetscMalloc(p*sizeof(int),&flz);
1776:     PetscMemcpy(flz,lz,p*sizeof(int));
1777:   }
1778:   da->lx = flx;
1779:   da->ly = fly;
1780:   da->lz = flz;

1782:   /*
1783:      Note the following will be removed soon. Since the functionality 
1784:     is replaced by the above. */

1786:   /* Construct the mapping from current global ordering to global
1787:      ordering that would be used if only 1 processor were employed.
1788:      This mapping is intended only for internal use by discrete
1789:      function and matrix viewers.

1791:      Note: At this point, x has already been adjusted for multiple
1792:      degrees of freedom per node.
1793:    */
1794:   ldim = x*y*z;
1795:   VecGetSize(global,&gdim);
1796:   PetscMalloc(gdim*sizeof(int),&da->gtog1);
1797:   PetscLogObjectMemory(da,gdim*sizeof(int));
1798:   PetscMalloc((2*(gdim+ldim))*sizeof(int),&gA);
1799:   gB        = (int *)(gA + ldim);
1800:   gAall     = (int *)(gB + ldim);
1801:   gBall     = (int *)(gAall + gdim);
1802:   /* Compute local parts of global orderings */
1803:   ict = 0;
1804:   for (k=zs; k<ze; k++) {
1805:     for (j=ys; j<ye; j++) {
1806:       for (i=xs; i<xe; i++) {
1807:         /* gA = global number for 1 proc; gB = current global number */
1808:         gA[ict] = i + j*M*dof + k*M*N*dof;
1809:         gB[ict] = start + ict;
1810:         ict++;
1811:       }
1812:     }
1813:   }
1814:   /* Broadcast the orderings */
1815:   MPI_Allgatherv(gA,ldim,MPI_INT,gAall,ldims,bases,MPI_INT,comm);
1816:   MPI_Allgatherv(gB,ldim,MPI_INT,gBall,ldims,bases,MPI_INT,comm);
1817:   for (i=0; i<gdim; i++) da->gtog1[gBall[i]] = gAall[i];
1818:   PetscFree(gA);
1819:   PetscFree(bases);

1821:   PetscOptionsHasName(PETSC_NULL,"-da_view",&flg1);
1822:   if (flg1) {DAView(da,PETSC_VIEWER_STDOUT_(da->comm));}
1823:   PetscOptionsHasName(PETSC_NULL,"-da_view_draw",&flg1);
1824:   if (flg1) {DAView(da,PETSC_VIEWER_DRAW_(da->comm));}
1825:   PetscOptionsHasName(PETSC_NULL,"-help",&flg1);
1826:   if (flg1) {DAPrintHelp(da);}
1827:   PetscPublishAll(da);

1829: #if defined(PETSC_HAVE_AMS)
1830:   PetscObjectComposeFunctionDynamic((PetscObject)global,"AMSSetFieldBlock_C",
1831:          "AMSSetFieldBlock_DA",AMSSetFieldBlock_DA);
1832:   PetscObjectComposeFunctionDynamic((PetscObject)local,"AMSSetFieldBlock_C",
1833:          "AMSSetFieldBlock_DA",AMSSetFieldBlock_DA);
1834:   if (((PetscObject)global)->amem > -1) {
1835:     AMSSetFieldBlock_DA(((PetscObject)global)->amem,"values",global);
1836:   }
1837: #endif
1838:   VecSetOperation(global,VECOP_VIEW,(void(*)())VecView_MPI_DA);
1839:   VecSetOperation(global,VECOP_LOADINTOVECTOR,(void(*)())VecLoadIntoVector_Binary_DA);
1840:   return(0);
1841: }