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
183: -  -da_noao - do not compute natural to PETSc ordering object

185:    Level: beginner

187:    Notes:
188:    If you are having problems with running out of memory than run with the option -da_noao

190:    The stencil type DA_STENCIL_STAR with width 1 corresponds to the 
191:    standard 7-pt stencil, while DA_STENCIL_BOX with width 1 denotes
192:    the standard 27-pt stencil.

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

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

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

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

226:   *inra = 0;
227: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
228:   DMInitializePackage(PETSC_NULL);
229: #endif

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

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

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

264:   PetscLogObjectCreate(da);
265:   PetscLogObjectMemory(da,sizeof(struct _p_DA));
266:   da->dim        = 3;
267:   da->interptype = DA_Q1;
268:   da->gtog1      = 0;
269:   da->refine_x   = refine_x;
270:   da->refine_y   = refine_y;
271:   da->refine_z   = refine_z;
272:   PetscMalloc(dof*sizeof(char*),&da->fieldname);
273:   PetscMemzero(da->fieldname,dof*sizeof(char*));

275:   MPI_Comm_size(comm,&size);
276:   MPI_Comm_rank(comm,&rank);

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

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

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

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

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

428:   /* X Periodic */
429:   if (DAXPeriodic(wrap)){
430:     Xs = xs - s;
431:     Xe = xe + s;
432:   }

434:   /* Y Periodic */
435:   if (DAYPeriodic(wrap)){
436:     Ys = ys - s;
437:     Ye = ye + s;
438:   }

440:   /* Z Periodic */
441:   if (DAZPeriodic(wrap)){
442:     Zs = zs - s;
443:     Ze = ze + s;
444:   }

446:   /* Resize all X parameters to reflect w */
447:   x   *= dof;
448:   xs  *= dof;
449:   xe  *= dof;
450:   Xs  *= dof;
451:   Xe  *= dof;
452:   s_x  = s*dof;
453:   s_y  = s;
454:   s_z  = s;

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

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

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

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

496:   VecScatterCreate(local,from,global,to,&ltog);
497:   PetscLogObjectParent(da,to);
498:   PetscLogObjectParent(da,from);
499:   PetscLogObjectParent(da,ltog);
500:   ISDestroy(from);
501:   ISDestroy(to);

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

547:   /* determine who lies on each side of use stored in    n24 n25 n26
548:                                                          n21 n22 n23
549:                                                          n18 n19 n20

551:                                                          n15 n16 n17
552:                                                          n12     n14
553:                                                          n9  n10 n11

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

574:   n9  = rank - m - 1;
575:   n10 = rank - m;
576:   n11 = rank - m + 1;
577:   n12 = rank - 1;
578:   n14 = rank + 1;
579:   n15 = rank + m - 1;
580:   n16 = rank + m;
581:   n17 = rank + m + 1;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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


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

787:   nn = 0;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1030:   PetscLogObjectParent(da,global);
1031:   PetscLogObjectParent(da,local);

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

1043:     nn = 0;

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

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

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

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

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

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

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

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

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

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

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

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

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

1247:       for (i=1; i<=s_y; i++) {
1248:         if (n24 >= 0) { /* left above */
1249:           x_t = lx[n24 % m]*dof;
1250:           y_t = ly[(n24 % (m*n))/m];
1251:           /* z_t = lz[n24 / (m*n)]; */
1252:           s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
1253:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1254:         }
1255:         if (n25 >= 0) { /* directly above */
1256:           x_t = x;
1257:           y_t = ly[(n25 % (m*n))/m];
1258:           /* z_t = lz[n25 / (m*n)]; */
1259:           s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
1260:           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1261:         }
1262:         if (n26 >= 0) { /* right above */
1263:           x_t = lx[n26 % m]*dof;
1264:           y_t = ly[(n26 % (m*n))/m];
1265:           /* z_t = lz[n26 / (m*n)]; */
1266:           s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
1267:           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1268:         }
1269:       }
1270:     }
1271:   }
1272:   da->global    = global;
1273:   da->local     = local;
1274:   da->gtol      = gtol;
1275:   da->ltog      = ltog;
1276:   da->idx       = idx;
1277:   da->Nl        = nn;
1278:   da->base      = base;
1279:   da->ops->view = DAView_3d;
1280:   da->wrap      = wrap;
1281:   *inra = da;

1283:   /* 
1284:      Set the local to global ordering in the global vector, this allows use
1285:      of VecSetValuesLocal().
1286:   */
1287:   ierr  = ISLocalToGlobalMappingCreate(comm,nn,idx,&da->ltogmap);
1288:   ierr  = VecSetLocalToGlobalMapping(da->global,da->ltogmap);
1289:   ISLocalToGlobalMappingBlock(da->ltogmap,da->w,&da->ltogmapb);
1290:   VecSetLocalToGlobalMappingBlock(da->global,da->ltogmapb);
1291:   PetscLogObjectParent(da,da->ltogmap);

1293:   /* redo idx to include "missing" ghost points */
1294:   /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
1295: 
1296:   /* Assume Nodes are Internal to the Cube */
1297: 
1298:   n0  = rank - m*n - m - 1;
1299:   n1  = rank - m*n - m;
1300:   n2  = rank - m*n - m + 1;
1301:   n3  = rank - m*n -1;
1302:   n4  = rank - m*n;
1303:   n5  = rank - m*n + 1;
1304:   n6  = rank - m*n + m - 1;
1305:   n7  = rank - m*n + m;
1306:   n8  = rank - m*n + m + 1;

1308:   n9  = rank - m - 1;
1309:   n10 = rank - m;
1310:   n11 = rank - m + 1;
1311:   n12 = rank - 1;
1312:   n14 = rank + 1;
1313:   n15 = rank + m - 1;
1314:   n16 = rank + m;
1315:   n17 = rank + m + 1;

1317:   n18 = rank + m*n - m - 1;
1318:   n19 = rank + m*n - m;
1319:   n20 = rank + m*n - m + 1;
1320:   n21 = rank + m*n - 1;
1321:   n22 = rank + m*n;
1322:   n23 = rank + m*n + 1;
1323:   n24 = rank + m*n + m - 1;
1324:   n25 = rank + m*n + m;
1325:   n26 = rank + m*n + m + 1;

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

1329:   if (xs == 0) { /* First assume not corner or edge */
1330:     n0  = rank       -1 - (m*n);
1331:     n3  = rank + m   -1 - (m*n);
1332:     n6  = rank + 2*m -1 - (m*n);
1333:     n9  = rank       -1;
1334:     n12 = rank + m   -1;
1335:     n15 = rank + 2*m -1;
1336:     n18 = rank       -1 + (m*n);
1337:     n21 = rank + m   -1 + (m*n);
1338:     n24 = rank + 2*m -1 + (m*n);
1339:    }

1341:   if (xe == M*dof) { /* First assume not corner or edge */
1342:     n2  = rank -2*m +1 - (m*n);
1343:     n5  = rank - m  +1 - (m*n);
1344:     n8  = rank      +1 - (m*n);
1345:     n11 = rank -2*m +1;
1346:     n14 = rank - m  +1;
1347:     n17 = rank      +1;
1348:     n20 = rank -2*m +1 + (m*n);
1349:     n23 = rank - m  +1 + (m*n);
1350:     n26 = rank      +1 + (m*n);
1351:   }

1353:   if (ys==0) { /* First assume not corner or edge */
1354:     n0  = rank + m * (n-1) -1 - (m*n);
1355:     n1  = rank + m * (n-1)    - (m*n);
1356:     n2  = rank + m * (n-1) +1 - (m*n);
1357:     n9  = rank + m * (n-1) -1;
1358:     n10 = rank + m * (n-1);
1359:     n11 = rank + m * (n-1) +1;
1360:     n18 = rank + m * (n-1) -1 + (m*n);
1361:     n19 = rank + m * (n-1)    + (m*n);
1362:     n20 = rank + m * (n-1) +1 + (m*n);
1363:   }

1365:   if (ye == N) { /* First assume not corner or edge */
1366:     n6  = rank - m * (n-1) -1 - (m*n);
1367:     n7  = rank - m * (n-1)    - (m*n);
1368:     n8  = rank - m * (n-1) +1 - (m*n);
1369:     n15 = rank - m * (n-1) -1;
1370:     n16 = rank - m * (n-1);
1371:     n17 = rank - m * (n-1) +1;
1372:     n24 = rank - m * (n-1) -1 + (m*n);
1373:     n25 = rank - m * (n-1)    + (m*n);
1374:     n26 = rank - m * (n-1) +1 + (m*n);
1375:   }
1376: 
1377:   if (zs == 0) { /* First assume not corner or edge */
1378:     n0 = size - (m*n) + rank - m - 1;
1379:     n1 = size - (m*n) + rank - m;
1380:     n2 = size - (m*n) + rank - m + 1;
1381:     n3 = size - (m*n) + rank - 1;
1382:     n4 = size - (m*n) + rank;
1383:     n5 = size - (m*n) + rank + 1;
1384:     n6 = size - (m*n) + rank + m - 1;
1385:     n7 = size - (m*n) + rank + m ;
1386:     n8 = size - (m*n) + rank + m + 1;
1387:   }

1389:   if (ze == P) { /* First assume not corner or edge */
1390:     n18 = (m*n) - (size-rank) - m - 1;
1391:     n19 = (m*n) - (size-rank) - m;
1392:     n20 = (m*n) - (size-rank) - m + 1;
1393:     n21 = (m*n) - (size-rank) - 1;
1394:     n22 = (m*n) - (size-rank);
1395:     n23 = (m*n) - (size-rank) + 1;
1396:     n24 = (m*n) - (size-rank) + m - 1;
1397:     n25 = (m*n) - (size-rank) + m;
1398:     n26 = (m*n) - (size-rank) + m + 1;
1399:   }

1401:   if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */
1402:     n0 = size - m*n + rank + m-1 - m;
1403:     n3 = size - m*n + rank + m-1;
1404:     n6 = size - m*n + rank + m-1 + m;
1405:   }
1406: 
1407:   if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */
1408:     n18 = m*n - (size - rank) + m-1 - m;
1409:     n21 = m*n - (size - rank) + m-1;
1410:     n24 = m*n - (size - rank) + m-1 + m;
1411:   }

1413:   if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */
1414:     n0  = rank + m*n -1 - m*n;
1415:     n9  = rank + m*n -1;
1416:     n18 = rank + m*n -1 + m*n;
1417:   }

1419:   if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */
1420:     n6  = rank - m*(n-1) + m-1 - m*n;
1421:     n15 = rank - m*(n-1) + m-1;
1422:     n24 = rank - m*(n-1) + m-1 + m*n;
1423:   }

1425:   if ((xe==M*dof) && (zs==0)) { /* Assume an edge, not corner */
1426:     n2 = size - (m*n-rank) - (m-1) - m;
1427:     n5 = size - (m*n-rank) - (m-1);
1428:     n8 = size - (m*n-rank) - (m-1) + m;
1429:   }

1431:   if ((xe==M*dof) && (ze==P)) { /* Assume an edge, not corner */
1432:     n20 = m*n - (size - rank) - (m-1) - m;
1433:     n23 = m*n - (size - rank) - (m-1);
1434:     n26 = m*n - (size - rank) - (m-1) + m;
1435:   }

1437:   if ((xe==M*dof) && (ys==0)) { /* Assume an edge, not corner */
1438:     n2  = rank + m*(n-1) - (m-1) - m*n;
1439:     n11 = rank + m*(n-1) - (m-1);
1440:     n20 = rank + m*(n-1) - (m-1) + m*n;
1441:   }

1443:   if ((xe==M*dof) && (ye==N)) { /* Assume an edge, not corner */
1444:     n8  = rank - m*n +1 - m*n;
1445:     n17 = rank - m*n +1;
1446:     n26 = rank - m*n +1 + m*n;
1447:   }

1449:   if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */
1450:     n0 = size - m + rank -1;
1451:     n1 = size - m + rank;
1452:     n2 = size - m + rank +1;
1453:   }

1455:   if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */
1456:     n18 = m*n - (size - rank) + m*(n-1) -1;
1457:     n19 = m*n - (size - rank) + m*(n-1);
1458:     n20 = m*n - (size - rank) + m*(n-1) +1;
1459:   }

1461:   if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */
1462:     n6 = size - (m*n-rank) - m * (n-1) -1;
1463:     n7 = size - (m*n-rank) - m * (n-1);
1464:     n8 = size - (m*n-rank) - m * (n-1) +1;
1465:   }

1467:   if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */
1468:     n24 = rank - (size-m) -1;
1469:     n25 = rank - (size-m);
1470:     n26 = rank - (size-m) +1;
1471:   }

1473:   /* Check for Corners */
1474:   if ((xs==0)   && (ys==0) && (zs==0)) { n0  = size -1;}
1475:   if ((xs==0)   && (ys==0) && (ze==P)) { n18 = m*n-1;}
1476:   if ((xs==0)   && (ye==N) && (zs==0)) { n6  = (size-1)-m*(n-1);}
1477:   if ((xs==0)   && (ye==N) && (ze==P)) { n24 = m-1;}
1478:   if ((xe==M*dof) && (ys==0) && (zs==0)) { n2  = size-m;}
1479:   if ((xe==M*dof) && (ys==0) && (ze==P)) { n20 = m*n-m;}
1480:   if ((xe==M*dof) && (ye==N) && (zs==0)) { n8  = size-m*n;}
1481:   if ((xe==M*dof) && (ye==N) && (ze==P)) { n26 = 0;}

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

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

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

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

1503:   nn = 0;

1505:   /* Bottom Level */
1506:   for (k=0; k<s_z; k++) {
1507:     for (i=1; i<=s_y; i++) {
1508:       if (n0 >= 0) { /* left below */
1509:         x_t = lx[n0 % m]*dof;
1510:         y_t = ly[(n0 % (m*n))/m];
1511:         z_t = lz[n0 / (m*n)];
1512:         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;
1513:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1514:       }
1515:       if (n1 >= 0) { /* directly below */
1516:         x_t = x;
1517:         y_t = ly[(n1 % (m*n))/m];
1518:         z_t = lz[n1 / (m*n)];
1519:         s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
1520:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1521:       }
1522:       if (n2 >= 0) { /* right below */
1523:         x_t = lx[n2 % m]*dof;
1524:         y_t = ly[(n2 % (m*n))/m];
1525:         z_t = lz[n2 / (m*n)];
1526:         s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
1527:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1528:       }
1529:     }

1531:     for (i=0; i<y; i++) {
1532:       if (n3 >= 0) { /* directly left */
1533:         x_t = lx[n3 % m]*dof;
1534:         y_t = y;
1535:         z_t = lz[n3 / (m*n)];
1536:         s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1537:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1538:       }

1540:       if (n4 >= 0) { /* middle */
1541:         x_t = x;
1542:         y_t = y;
1543:         z_t = lz[n4 / (m*n)];
1544:         s_t = bases[n4] + i*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:       }

1548:       if (n5 >= 0) { /* directly right */
1549:         x_t = lx[n5 % m]*dof;
1550:         y_t = y;
1551:         z_t = lz[n5 / (m*n)];
1552:         s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1553:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1554:       }
1555:     }

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

1582:   /* Middle Level */
1583:   for (k=0; k<z; k++) {
1584:     for (i=1; i<=s_y; i++) {
1585:       if (n9 >= 0) { /* left below */
1586:         x_t = lx[n9 % m]*dof;
1587:         y_t = ly[(n9 % (m*n))/m];
1588:         /* z_t = z; */
1589:         s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
1590:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1591:       }
1592:       if (n10 >= 0) { /* directly below */
1593:         x_t = x;
1594:         y_t = ly[(n10 % (m*n))/m];
1595:         /* z_t = z; */
1596:         s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1597:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1598:       }
1599:       if (n11 >= 0) { /* right below */
1600:         x_t = lx[n11 % m]*dof;
1601:         y_t = ly[(n11 % (m*n))/m];
1602:         /* z_t = z; */
1603:         s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1604:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1605:       }
1606:     }

1608:     for (i=0; i<y; i++) {
1609:       if (n12 >= 0) { /* directly left */
1610:         x_t = lx[n12 % m]*dof;
1611:         y_t = y;
1612:         /* z_t = z; */
1613:         s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
1614:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1615:       }

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

1621:       if (n14 >= 0) { /* directly right */
1622:         x_t = lx[n14 % m]*dof;
1623:         y_t = y;
1624:         /* z_t = z; */
1625:         s_t = bases[n14] + i*x_t + k*x_t*y_t;
1626:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1627:       }
1628:     }

1630:     for (i=1; i<=s_y; i++) {
1631:       if (n15 >= 0) { /* left above */
1632:         x_t = lx[n15 % m]*dof;
1633:         y_t = ly[(n15 % (m*n))/m];
1634:         /* z_t = z; */
1635:         s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
1636:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1637:       }
1638:       if (n16 >= 0) { /* directly above */
1639:         x_t = x;
1640:         y_t = ly[(n16 % (m*n))/m];
1641:         /* z_t = z; */
1642:         s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
1643:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1644:       }
1645:       if (n17 >= 0) { /* right above */
1646:         x_t = lx[n17 % m]*dof;
1647:         y_t = ly[(n17 % (m*n))/m];
1648:         /* z_t = z; */
1649:         s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
1650:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1651:       }
1652:     }
1653:   }
1654: 
1655:   /* Upper Level */
1656:   for (k=0; k<s_z; k++) {
1657:     for (i=1; i<=s_y; i++) {
1658:       if (n18 >= 0) { /* left below */
1659:         x_t = lx[n18 % m]*dof;
1660:         y_t = ly[(n18 % (m*n))/m];
1661:         /* z_t = lz[n18 / (m*n)]; */
1662:         s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
1663:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1664:       }
1665:       if (n19 >= 0) { /* directly below */
1666:         x_t = x;
1667:         y_t = ly[(n19 % (m*n))/m];
1668:         /* z_t = lz[n19 / (m*n)]; */
1669:         s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1670:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1671:       }
1672:       if (n20 >= 0) { /* right belodof */
1673:         x_t = lx[n20 % m]*dof;
1674:         y_t = ly[(n20 % (m*n))/m];
1675:         /* z_t = lz[n20 / (m*n)]; */
1676:         s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1677:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1678:       }
1679:     }

1681:     for (i=0; i<y; i++) {
1682:       if (n21 >= 0) { /* directly left */
1683:         x_t = lx[n21 % m]*dof;
1684:         y_t = y;
1685:         /* z_t = lz[n21 / (m*n)]; */
1686:         s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
1687:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1688:       }

1690:       if (n22 >= 0) { /* middle */
1691:         x_t = x;
1692:         y_t = y;
1693:         /* z_t = lz[n22 / (m*n)]; */
1694:         s_t = bases[n22] + i*x_t + k*x_t*y_t;
1695:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1696:       }

1698:       if (n23 >= 0) { /* directly right */
1699:         x_t = lx[n23 % m]*dof;
1700:         y_t = y;
1701:         /* z_t = lz[n23 / (m*n)]; */
1702:         s_t = bases[n23] + i*x_t + k*x_t*y_t;
1703:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1704:       }
1705:     }

1707:     for (i=1; i<=s_y; i++) {
1708:       if (n24 >= 0) { /* left above */
1709:         x_t = lx[n24 % m]*dof;
1710:         y_t = ly[(n24 % (m*n))/m];
1711:         /* z_t = lz[n24 / (m*n)]; */
1712:         s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
1713:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1714:       }
1715:       if (n25 >= 0) { /* directly above */
1716:         x_t = x;
1717:         y_t = ly[(n25 % (m*n))/m];
1718:         /* z_t = lz[n25 / (m*n)]; */
1719:         s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
1720:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1721:       }
1722:       if (n26 >= 0) { /* right above */
1723:         x_t = lx[n26 % m]*dof;
1724:         y_t = ly[(n26 % (m*n))/m];
1725:         /* z_t = lz[n26 / (m*n)]; */
1726:         s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
1727:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1728:       }
1729:     }
1730:   }

1732:   /* construct the local to local scatter context */
1733:   /* 
1734:       We simply remap the values in the from part of 
1735:     global to local to read from an array with the ghost values 
1736:     rather then from the plan array.
1737:   */
1738:   VecScatterCopy(gtol,&da->ltol);
1739:   PetscLogObjectParent(da,da->ltol);
1740:   left   = xs - Xs;
1741:   bottom = ys - Ys; top = bottom + y;
1742:   down   = zs - Zs; up  = down + z;
1743:   count  = x*(top-bottom)*(up-down);
1744:   PetscMalloc(count*sizeof(int),&idx);
1745:   count  = 0;
1746:   for (i=down; i<up; i++) {
1747:     for (j=bottom; j<top; j++) {
1748:       for (k=0; k<x; k++) {
1749:         idx[count++] = (left+j*(Xe-Xs))+i*(Xe-Xs)*(Ye-Ys) + k;
1750:       }
1751:     }
1752:   }
1753:   VecScatterRemap(da->ltol,idx,PETSC_NULL);
1754:   PetscFree(idx);


1757:   /* 
1758:      Build the natural ordering to PETSc ordering mappings.
1759:   */
1760:   PetscOptionsHasName(PETSC_NULL,"-da_noao",&flg1);
1761:   if (!flg1) {
1762:     IS  ispetsc,isnatural;
1763:     int *lidx,lict = 0;
1764:     int Nlocal = (da->xe-da->xs)*(da->ye-da->ys)*(da->ze-da->zs);

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

1768:     PetscMalloc(Nlocal*sizeof(int),&lidx);
1769:     for (k=zs; k<ze; k++) {
1770:       for (j=ys; j<ye; j++) {
1771:         for (i=xs; i<xe; i++) {
1772:           lidx[lict++] = i + j*M*dof + k*M*N*dof;
1773:         }
1774:       }
1775:     }
1776:     ISCreateGeneral(comm,Nlocal,lidx,&isnatural);
1777:     PetscFree(lidx);

1779:     AOCreateBasicIS(isnatural,ispetsc,&da->ao);
1780:     PetscLogObjectParent(da,da->ao);
1781:     ISDestroy(ispetsc);
1782:     ISDestroy(isnatural);
1783:   } else {
1784:     da->ao = PETSC_NULL;
1785:   }

1787:   if (!flx) {
1788:     PetscMalloc(m*sizeof(int),&flx);
1789:     PetscMemcpy(flx,lx,m*sizeof(int));
1790:   }
1791:   if (!fly) {
1792:     PetscMalloc(n*sizeof(int),&fly);
1793:     PetscMemcpy(fly,ly,n*sizeof(int));
1794:   }
1795:   if (!flz) {
1796:     PetscMalloc(p*sizeof(int),&flz);
1797:     PetscMemcpy(flz,lz,p*sizeof(int));
1798:   }
1799:   da->lx = flx;
1800:   da->ly = fly;
1801:   da->lz = flz;

1803:   /*
1804:      Note the following will be removed soon. Since the functionality 
1805:     is replaced by the above. */

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

1812:      Note: At this point, x has already been adjusted for multiple
1813:      degrees of freedom per node.
1814:    */
1815:   ldim = x*y*z;
1816:   VecGetSize(global,&gdim);
1817:   PetscMalloc(gdim*sizeof(int),&da->gtog1);
1818:   PetscLogObjectMemory(da,gdim*sizeof(int));
1819:   PetscMalloc((2*(gdim+ldim))*sizeof(int),&gA);
1820:   gB        = (int *)(gA + ldim);
1821:   gAall     = (int *)(gB + ldim);
1822:   gBall     = (int *)(gAall + gdim);
1823:   /* Compute local parts of global orderings */
1824:   ict = 0;
1825:   for (k=zs; k<ze; k++) {
1826:     for (j=ys; j<ye; j++) {
1827:       for (i=xs; i<xe; i++) {
1828:         /* gA = global number for 1 proc; gB = current global number */
1829:         gA[ict] = i + j*M*dof + k*M*N*dof;
1830:         gB[ict] = start + ict;
1831:         ict++;
1832:       }
1833:     }
1834:   }
1835:   /* Broadcast the orderings */
1836:   MPI_Allgatherv(gA,ldim,MPI_INT,gAall,ldims,bases,MPI_INT,comm);
1837:   MPI_Allgatherv(gB,ldim,MPI_INT,gBall,ldims,bases,MPI_INT,comm);
1838:   for (i=0; i<gdim; i++) da->gtog1[gBall[i]] = gAall[i];
1839:   PetscFree(gA);
1840:   PetscFree(bases);

1842:   PetscOptionsHasName(PETSC_NULL,"-da_view",&flg1);
1843:   if (flg1) {DAView(da,PETSC_VIEWER_STDOUT_(da->comm));}
1844:   PetscOptionsHasName(PETSC_NULL,"-da_view_draw",&flg1);
1845:   if (flg1) {DAView(da,PETSC_VIEWER_DRAW_(da->comm));}
1846:   PetscOptionsHasName(PETSC_NULL,"-help",&flg1);
1847:   if (flg1) {DAPrintHelp(da);}
1848:   PetscPublishAll(da);

1850: #if defined(PETSC_HAVE_AMS)
1851:   PetscObjectComposeFunctionDynamic((PetscObject)global,"AMSSetFieldBlock_C",
1852:          "AMSSetFieldBlock_DA",AMSSetFieldBlock_DA);
1853:   PetscObjectComposeFunctionDynamic((PetscObject)local,"AMSSetFieldBlock_C",
1854:          "AMSSetFieldBlock_DA",AMSSetFieldBlock_DA);
1855:   if (((PetscObject)global)->amem > -1) {
1856:     AMSSetFieldBlock_DA(((PetscObject)global)->amem,"values",global);
1857:   }
1858: #endif
1859:   VecSetOperation(global,VECOP_VIEW,(void(*)(void))VecView_MPI_DA);
1860:   VecSetOperation(global,VECOP_LOADINTOVECTOR,(void(*)(void))VecLoadIntoVector_Binary_DA);
1861:   return(0);
1862: }