xref: /petsc/src/dm/impls/da/da3.c (revision 564755cd2444013b9df4620f3ca4ef575626a108)
147c6ae99SBarry Smith #define PETSCDM_DLL
247c6ae99SBarry Smith /*
347c6ae99SBarry Smith    Code for manipulating distributed regular 3d arrays in parallel.
447c6ae99SBarry Smith    File created by Peter Mell  7/14/95
547c6ae99SBarry Smith  */
647c6ae99SBarry Smith 
747c6ae99SBarry Smith #include "private/daimpl.h"     /*I   "petscda.h"    I*/
847c6ae99SBarry Smith 
947c6ae99SBarry Smith #undef __FUNCT__
109a42bb27SBarry Smith #define __FUNCT__ "DMView_DA_3d"
119a42bb27SBarry Smith PetscErrorCode DMView_DA_3d(DM da,PetscViewer viewer)
1247c6ae99SBarry Smith {
1347c6ae99SBarry Smith   PetscErrorCode ierr;
1447c6ae99SBarry Smith   PetscMPIInt    rank;
159a42bb27SBarry Smith   PetscBool      iascii,isdraw,isbinary;
1647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
179a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
189a42bb27SBarry Smith   PetscBool      ismatlab;
199a42bb27SBarry Smith #endif
2047c6ae99SBarry Smith 
2147c6ae99SBarry Smith   PetscFunctionBegin;
2247c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)da)->comm,&rank);CHKERRQ(ierr);
2347c6ae99SBarry Smith 
2447c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
2547c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
269a42bb27SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
279a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
289a42bb27SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr);
299a42bb27SBarry Smith #endif
3047c6ae99SBarry Smith   if (iascii) {
3147c6ae99SBarry Smith     PetscViewerFormat format;
3247c6ae99SBarry Smith 
3347c6ae99SBarry Smith     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
3447c6ae99SBarry Smith     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) {
3547c6ae99SBarry Smith       DALocalInfo info;
3647c6ae99SBarry Smith       ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
3747c6ae99SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %D N %D P %D m %D n %D p %D w %D s %D\n",rank,dd->M,dd->N,dd->P,dd->m,dd->n,dd->p,dd->w,dd->s);CHKERRQ(ierr);
3847c6ae99SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D, Y range of indices: %D %D, Z range of indices: %D %D\n",
3947c6ae99SBarry Smith                                                 info.xs,info.xs+info.xm,info.ys,info.ys+info.ym,info.zs,info.zs+info.zm);CHKERRQ(ierr);
4047c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX)
4147c6ae99SBarry Smith       if (dd->coordinates) {
4247c6ae99SBarry Smith         PetscInt        last;
4347c6ae99SBarry Smith         const PetscReal *coors;
4447c6ae99SBarry Smith         ierr = VecGetArrayRead(dd->coordinates,&coors);CHKERRQ(ierr);
4547c6ae99SBarry Smith         ierr = VecGetLocalSize(dd->coordinates,&last);CHKERRQ(ierr);
4647c6ae99SBarry Smith         last = last - 3;
4747c6ae99SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Lower left corner %G %G %G : Upper right %G %G %G\n",coors[0],coors[1],coors[2],coors[last],coors[last+1],coors[last+2]);CHKERRQ(ierr);
4847c6ae99SBarry Smith         ierr = VecRestoreArrayRead(dd->coordinates,&coors);CHKERRQ(ierr);
4947c6ae99SBarry Smith       }
5047c6ae99SBarry Smith #endif
5147c6ae99SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5247c6ae99SBarry Smith     }
5347c6ae99SBarry Smith   } else if (isdraw) {
5447c6ae99SBarry Smith     PetscDraw       draw;
5547c6ae99SBarry Smith     PetscReal     ymin = -1.0,ymax = (PetscReal)dd->N;
5647c6ae99SBarry Smith     PetscReal     xmin = -1.0,xmax = (PetscReal)((dd->M+2)*dd->P),x,y,ycoord,xcoord;
5747c6ae99SBarry Smith     PetscInt        k,plane,base,*idx;
5847c6ae99SBarry Smith     char       node[10];
5947c6ae99SBarry Smith     PetscBool  isnull;
6047c6ae99SBarry Smith 
6147c6ae99SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
6247c6ae99SBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
6347c6ae99SBarry Smith     ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
6447c6ae99SBarry Smith     ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
6547c6ae99SBarry Smith 
6647c6ae99SBarry Smith     /* first processor draw all node lines */
6747c6ae99SBarry Smith     if (!rank) {
6847c6ae99SBarry Smith       for (k=0; k<dd->P; k++) {
6947c6ae99SBarry Smith         ymin = 0.0; ymax = (PetscReal)(dd->N - 1);
7047c6ae99SBarry Smith         for (xmin=(PetscReal)(k*(dd->M+1)); xmin<(PetscReal)(dd->M+(k*(dd->M+1))); xmin++) {
7147c6ae99SBarry Smith           ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
7247c6ae99SBarry Smith         }
7347c6ae99SBarry Smith 
7447c6ae99SBarry Smith         xmin = (PetscReal)(k*(dd->M+1)); xmax = xmin + (PetscReal)(dd->M - 1);
7547c6ae99SBarry Smith         for (ymin=0; ymin<(PetscReal)dd->N; ymin++) {
7647c6ae99SBarry Smith           ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
7747c6ae99SBarry Smith         }
7847c6ae99SBarry Smith       }
7947c6ae99SBarry Smith     }
8047c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
8147c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
8247c6ae99SBarry Smith 
8347c6ae99SBarry Smith     for (k=0; k<dd->P; k++) {  /*Go through and draw for each plane*/
8447c6ae99SBarry Smith       if ((k >= dd->zs) && (k < dd->ze)) {
8547c6ae99SBarry Smith         /* draw my box */
8647c6ae99SBarry Smith         ymin = dd->ys;
8747c6ae99SBarry Smith         ymax = dd->ye - 1;
8847c6ae99SBarry Smith         xmin = dd->xs/dd->w    + (dd->M+1)*k;
8947c6ae99SBarry Smith         xmax =(dd->xe-1)/dd->w + (dd->M+1)*k;
9047c6ae99SBarry Smith 
9147c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
9247c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
9347c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
9447c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
9547c6ae99SBarry Smith 
9647c6ae99SBarry Smith         xmin = dd->xs/dd->w;
9747c6ae99SBarry Smith         xmax =(dd->xe-1)/dd->w;
9847c6ae99SBarry Smith 
9947c6ae99SBarry Smith         /* put in numbers*/
10047c6ae99SBarry Smith         base = (dd->base+(dd->xe-dd->xs)*(dd->ye-dd->ys)*(k-dd->zs))/dd->w;
10147c6ae99SBarry Smith 
10247c6ae99SBarry Smith         /* Identify which processor owns the box */
10347c6ae99SBarry Smith         sprintf(node,"%d",rank);
10447c6ae99SBarry Smith         ierr = PetscDrawString(draw,xmin+(dd->M+1)*k+.2,ymin+.3,PETSC_DRAW_RED,node);CHKERRQ(ierr);
10547c6ae99SBarry Smith 
10647c6ae99SBarry Smith         for (y=ymin; y<=ymax; y++) {
10747c6ae99SBarry Smith           for (x=xmin+(dd->M+1)*k; x<=xmax+(dd->M+1)*k; x++) {
10847c6ae99SBarry Smith             sprintf(node,"%d",(int)base++);
10947c6ae99SBarry Smith             ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
11047c6ae99SBarry Smith           }
11147c6ae99SBarry Smith         }
11247c6ae99SBarry Smith 
11347c6ae99SBarry Smith       }
11447c6ae99SBarry Smith     }
11547c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
11647c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
11747c6ae99SBarry Smith 
11847c6ae99SBarry Smith     for (k=0-dd->s; k<dd->P+dd->s; k++) {
11947c6ae99SBarry Smith       /* Go through and draw for each plane */
12047c6ae99SBarry Smith       if ((k >= dd->Zs) && (k < dd->Ze)) {
12147c6ae99SBarry Smith 
12247c6ae99SBarry Smith         /* overlay ghost numbers, useful for error checking */
12347c6ae99SBarry Smith         base = (dd->Xe-dd->Xs)*(dd->Ye-dd->Ys)*(k-dd->Zs); idx = dd->idx;
12447c6ae99SBarry Smith         plane=k;
12547c6ae99SBarry Smith         /* Keep z wrap around points on the dradrawg */
12647c6ae99SBarry Smith         if (k<0)    { plane=dd->P+k; }
12747c6ae99SBarry Smith         if (k>=dd->P) { plane=k-dd->P; }
12847c6ae99SBarry Smith         ymin = dd->Ys; ymax = dd->Ye;
12947c6ae99SBarry Smith         xmin = (dd->M+1)*plane*dd->w;
13047c6ae99SBarry Smith         xmax = (dd->M+1)*plane*dd->w+dd->M*dd->w;
13147c6ae99SBarry Smith         for (y=ymin; y<ymax; y++) {
13247c6ae99SBarry Smith           for (x=xmin+dd->Xs; x<xmin+dd->Xe; x+=dd->w) {
13347c6ae99SBarry Smith             sprintf(node,"%d",(int)(idx[base]/dd->w));
13447c6ae99SBarry Smith             ycoord = y;
13547c6ae99SBarry Smith             /*Keep y wrap around points on drawing */
13647c6ae99SBarry Smith             if (y<0)      { ycoord = dd->N+y; }
13747c6ae99SBarry Smith 
13847c6ae99SBarry Smith             if (y>=dd->N) { ycoord = y-dd->N; }
13947c6ae99SBarry Smith             xcoord = x;   /* Keep x wrap points on drawing */
14047c6ae99SBarry Smith 
14147c6ae99SBarry Smith             if (x<xmin)  { xcoord = xmax - (xmin-x); }
14247c6ae99SBarry Smith             if (x>=xmax) { xcoord = xmin + (x-xmax); }
14347c6ae99SBarry Smith             ierr = PetscDrawString(draw,xcoord/dd->w,ycoord,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
14447c6ae99SBarry Smith             base+=dd->w;
14547c6ae99SBarry Smith           }
14647c6ae99SBarry Smith         }
14747c6ae99SBarry Smith       }
14847c6ae99SBarry Smith     }
14947c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
15047c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
1519a42bb27SBarry Smith   } else if (isbinary){
1529a42bb27SBarry Smith     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
1539a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1549a42bb27SBarry Smith   } else if (ismatlab) {
1559a42bb27SBarry Smith     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
1569a42bb27SBarry Smith #endif
1579a42bb27SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for DA 1d",((PetscObject)viewer)->type_name);
15847c6ae99SBarry Smith   PetscFunctionReturn(0);
15947c6ae99SBarry Smith }
16047c6ae99SBarry Smith 
16147c6ae99SBarry Smith #undef __FUNCT__
1629a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_3D"
1639a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMSetUp_DA_3D(DM da)
16447c6ae99SBarry Smith {
16547c6ae99SBarry Smith   DM_DA               *dd           = (DM_DA*)da->data;
16647c6ae99SBarry Smith   const PetscInt       M            = dd->M;
16747c6ae99SBarry Smith   const PetscInt       N            = dd->N;
16847c6ae99SBarry Smith   const PetscInt       P            = dd->P;
16947c6ae99SBarry Smith   PetscInt             m            = dd->m;
17047c6ae99SBarry Smith   PetscInt             n            = dd->n;
17147c6ae99SBarry Smith   PetscInt             p            = dd->p;
17247c6ae99SBarry Smith   const PetscInt       dof          = dd->w;
17347c6ae99SBarry Smith   const PetscInt       s            = dd->s;
17447c6ae99SBarry Smith   const DAPeriodicType wrap         = dd->wrap;
17547c6ae99SBarry Smith   const DAStencilType  stencil_type = dd->stencil_type;
17647c6ae99SBarry Smith   PetscInt            *lx           = dd->lx;
17747c6ae99SBarry Smith   PetscInt            *ly           = dd->ly;
17847c6ae99SBarry Smith   PetscInt            *lz           = dd->lz;
17947c6ae99SBarry Smith   MPI_Comm             comm;
18047c6ae99SBarry Smith   PetscMPIInt    rank,size;
18147c6ae99SBarry Smith   PetscInt       xs = 0,xe,ys = 0,ye,zs = 0,ze,x = 0,y = 0,z = 0,Xs,Xe,Ys,Ye,Zs,Ze,start,end,pm;
18247c6ae99SBarry Smith   PetscInt       left,up,down,bottom,top,i,j,k,*idx,nn;
18347c6ae99SBarry Smith   PetscInt       n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14;
18447c6ae99SBarry Smith   PetscInt       n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26;
18547c6ae99SBarry Smith   PetscInt       *bases,*ldims,x_t,y_t,z_t,s_t,base,count,s_x,s_y,s_z;
18647c6ae99SBarry Smith   PetscInt       sn0 = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0;
18747c6ae99SBarry Smith   PetscInt       sn8 = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0;
18847c6ae99SBarry Smith   PetscInt       sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0;
18947c6ae99SBarry Smith   Vec            local,global;
19047c6ae99SBarry Smith   VecScatter     ltog,gtol;
19147c6ae99SBarry Smith   IS             to,from;
19247c6ae99SBarry Smith   PetscErrorCode ierr;
19347c6ae99SBarry Smith 
19447c6ae99SBarry Smith   PetscFunctionBegin;
19547c6ae99SBarry Smith   if (dof < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
19647c6ae99SBarry Smith   if (s < 0) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
19747c6ae99SBarry Smith 
19847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr);
19947c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
20047c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
20147c6ae99SBarry Smith 
20247c6ae99SBarry Smith   dd->dim = 3;
20347c6ae99SBarry Smith   ierr = PetscMalloc(dof*sizeof(char*),&dd->fieldname);CHKERRQ(ierr);
20447c6ae99SBarry Smith   ierr = PetscMemzero(dd->fieldname,dof*sizeof(char*));CHKERRQ(ierr);
20547c6ae99SBarry Smith 
20647c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
20747c6ae99SBarry Smith     if (m < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
20847c6ae99SBarry Smith     else if (m > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
20947c6ae99SBarry Smith   }
21047c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
21147c6ae99SBarry Smith     if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
21247c6ae99SBarry Smith     else if (n > size)  SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
21347c6ae99SBarry Smith   }
21447c6ae99SBarry Smith   if (p != PETSC_DECIDE) {
21547c6ae99SBarry Smith     if (p < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %D",p);
21647c6ae99SBarry Smith     else if (p > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %D %d",p,size);
21747c6ae99SBarry Smith   }
21847c6ae99SBarry Smith 
21947c6ae99SBarry Smith   /* Partition the array among the processors */
22047c6ae99SBarry Smith   if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) {
22147c6ae99SBarry Smith     m = size/(n*p);
22247c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
22347c6ae99SBarry Smith     n = size/(m*p);
22447c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
22547c6ae99SBarry Smith     p = size/(m*n);
22647c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
22747c6ae99SBarry Smith     /* try for squarish distribution */
22847c6ae99SBarry Smith     m = (int)(0.5 + sqrt(((PetscReal)M)*((PetscReal)size)/((PetscReal)N*p)));
22947c6ae99SBarry Smith     if (!m) m = 1;
23047c6ae99SBarry Smith     while (m > 0) {
23147c6ae99SBarry Smith       n = size/(m*p);
23247c6ae99SBarry Smith       if (m*n*p == size) break;
23347c6ae99SBarry Smith       m--;
23447c6ae99SBarry Smith     }
23547c6ae99SBarry Smith     if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %D",p);
23647c6ae99SBarry Smith     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
23747c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
23847c6ae99SBarry Smith     /* try for squarish distribution */
23947c6ae99SBarry Smith     m = (int)(0.5 + sqrt(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n)));
24047c6ae99SBarry Smith     if (!m) m = 1;
24147c6ae99SBarry Smith     while (m > 0) {
24247c6ae99SBarry Smith       p = size/(m*n);
24347c6ae99SBarry Smith       if (m*n*p == size) break;
24447c6ae99SBarry Smith       m--;
24547c6ae99SBarry Smith     }
24647c6ae99SBarry Smith     if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %D",n);
24747c6ae99SBarry Smith     if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;}
24847c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
24947c6ae99SBarry Smith     /* try for squarish distribution */
25047c6ae99SBarry Smith     n = (int)(0.5 + sqrt(((PetscReal)N)*((PetscReal)size)/((PetscReal)P*m)));
25147c6ae99SBarry Smith     if (!n) n = 1;
25247c6ae99SBarry Smith     while (n > 0) {
25347c6ae99SBarry Smith       p = size/(m*n);
25447c6ae99SBarry Smith       if (m*n*p == size) break;
25547c6ae99SBarry Smith       n--;
25647c6ae99SBarry Smith     }
25747c6ae99SBarry Smith     if (!n) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %D",n);
25847c6ae99SBarry Smith     if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;}
25947c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
26047c6ae99SBarry Smith     /* try for squarish distribution */
26147c6ae99SBarry Smith     n = (PetscInt)(0.5 + pow(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),(PetscReal)(1./3.)));
26247c6ae99SBarry Smith     if (!n) n = 1;
26347c6ae99SBarry Smith     while (n > 0) {
26447c6ae99SBarry Smith       pm = size/n;
26547c6ae99SBarry Smith       if (n*pm == size) break;
26647c6ae99SBarry Smith       n--;
26747c6ae99SBarry Smith     }
26847c6ae99SBarry Smith     if (!n) n = 1;
26947c6ae99SBarry Smith     m = (PetscInt)(0.5 + sqrt(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n)));
27047c6ae99SBarry Smith     if (!m) m = 1;
27147c6ae99SBarry Smith     while (m > 0) {
27247c6ae99SBarry Smith       p = size/(m*n);
27347c6ae99SBarry Smith       if (m*n*p == size) break;
27447c6ae99SBarry Smith       m--;
27547c6ae99SBarry Smith     }
27647c6ae99SBarry Smith     if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;}
27747c6ae99SBarry Smith   } else if (m*n*p != size) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
27847c6ae99SBarry Smith 
27947c6ae99SBarry Smith   if (m*n*p != size) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_PLIB,"Could not find good partition");
28047c6ae99SBarry Smith   if (M < m) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
28147c6ae99SBarry Smith   if (N < n) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
28247c6ae99SBarry Smith   if (P < p) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %D %D",P,p);
28347c6ae99SBarry Smith 
28447c6ae99SBarry Smith   /*
28547c6ae99SBarry Smith      Determine locally owned region
28647c6ae99SBarry Smith      [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes
28747c6ae99SBarry Smith   */
28847c6ae99SBarry Smith 
28947c6ae99SBarry Smith   if (!lx) {
29047c6ae99SBarry Smith     ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
29147c6ae99SBarry Smith     lx = dd->lx;
29247c6ae99SBarry Smith     for (i=0; i<m; i++) {
29347c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > (i % m));
29447c6ae99SBarry Smith     }
29547c6ae99SBarry Smith   }
29647c6ae99SBarry Smith   x  = lx[rank % m];
29747c6ae99SBarry Smith   xs = 0;
29847c6ae99SBarry Smith   for (i=0; i<(rank%m); i++) { xs += lx[i];}
29947c6ae99SBarry Smith   if (m > 1 && x < s) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column width is too thin for stencil! %D %D",x,s);
30047c6ae99SBarry Smith 
30147c6ae99SBarry Smith   if (!ly) {
30247c6ae99SBarry Smith     ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
30347c6ae99SBarry Smith     ly = dd->ly;
30447c6ae99SBarry Smith     for (i=0; i<n; i++) {
30547c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > (i % n));
30647c6ae99SBarry Smith     }
30747c6ae99SBarry Smith   }
30847c6ae99SBarry Smith   y  = ly[(rank % (m*n))/m];
30947c6ae99SBarry Smith   if (n > 1 && y < s) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row width is too thin for stencil! %D %D",y,s);
31047c6ae99SBarry Smith   ys = 0;
31147c6ae99SBarry Smith   for (i=0; i<(rank % (m*n))/m; i++) { ys += ly[i];}
31247c6ae99SBarry Smith 
31347c6ae99SBarry Smith   if (!lz) {
31447c6ae99SBarry Smith     ierr = PetscMalloc(p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr);
31547c6ae99SBarry Smith     lz = dd->lz;
31647c6ae99SBarry Smith     for (i=0; i<p; i++) {
31747c6ae99SBarry Smith       lz[i] = P/p + ((P % p) > (i % p));
31847c6ae99SBarry Smith     }
31947c6ae99SBarry Smith   }
32047c6ae99SBarry Smith   z  = lz[rank/(m*n)];
32147c6ae99SBarry Smith   if (p > 1 && z < s) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Plane width is too thin for stencil! %D %D",z,s);
32247c6ae99SBarry Smith   zs = 0;
32347c6ae99SBarry Smith   for (i=0; i<(rank/(m*n)); i++) { zs += lz[i];}
32447c6ae99SBarry Smith   ye = ys + y;
32547c6ae99SBarry Smith   xe = xs + x;
32647c6ae99SBarry Smith   ze = zs + z;
32747c6ae99SBarry Smith 
32847c6ae99SBarry Smith   /* determine ghost region */
32947c6ae99SBarry Smith   /* Assume No Periodicity */
33047c6ae99SBarry Smith   if (xs-s > 0) Xs = xs - s; else Xs = 0;
33147c6ae99SBarry Smith   if (ys-s > 0) Ys = ys - s; else Ys = 0;
33247c6ae99SBarry Smith   if (zs-s > 0) Zs = zs - s; else Zs = 0;
33347c6ae99SBarry Smith   if (xe+s <= M) Xe = xe + s; else Xe = M;
33447c6ae99SBarry Smith   if (ye+s <= N) Ye = ye + s; else Ye = N;
33547c6ae99SBarry Smith   if (ze+s <= P) Ze = ze + s; else Ze = P;
33647c6ae99SBarry Smith 
33747c6ae99SBarry Smith   /* X Periodic */
33847c6ae99SBarry Smith   if (DAXPeriodic(wrap)){
33947c6ae99SBarry Smith     Xs = xs - s;
34047c6ae99SBarry Smith     Xe = xe + s;
34147c6ae99SBarry Smith   }
34247c6ae99SBarry Smith 
34347c6ae99SBarry Smith   /* Y Periodic */
34447c6ae99SBarry Smith   if (DAYPeriodic(wrap)){
34547c6ae99SBarry Smith     Ys = ys - s;
34647c6ae99SBarry Smith     Ye = ye + s;
34747c6ae99SBarry Smith   }
34847c6ae99SBarry Smith 
34947c6ae99SBarry Smith   /* Z Periodic */
35047c6ae99SBarry Smith   if (DAZPeriodic(wrap)){
35147c6ae99SBarry Smith     Zs = zs - s;
35247c6ae99SBarry Smith     Ze = ze + s;
35347c6ae99SBarry Smith   }
35447c6ae99SBarry Smith 
35547c6ae99SBarry Smith   /* Resize all X parameters to reflect w */
35647c6ae99SBarry Smith   x   *= dof;
35747c6ae99SBarry Smith   xs  *= dof;
35847c6ae99SBarry Smith   xe  *= dof;
35947c6ae99SBarry Smith   Xs  *= dof;
36047c6ae99SBarry Smith   Xe  *= dof;
36147c6ae99SBarry Smith   s_x  = s*dof;
36247c6ae99SBarry Smith   s_y  = s;
36347c6ae99SBarry Smith   s_z  = s;
36447c6ae99SBarry Smith 
36547c6ae99SBarry Smith   /* determine starting point of each processor */
36647c6ae99SBarry Smith   nn       = x*y*z;
36747c6ae99SBarry Smith   ierr     = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr);
36847c6ae99SBarry Smith   ierr     = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
36947c6ae99SBarry Smith   bases[0] = 0;
37047c6ae99SBarry Smith   for (i=1; i<=size; i++) {
37147c6ae99SBarry Smith     bases[i] = ldims[i-1];
37247c6ae99SBarry Smith   }
37347c6ae99SBarry Smith   for (i=1; i<=size; i++) {
37447c6ae99SBarry Smith     bases[i] += bases[i-1];
37547c6ae99SBarry Smith   }
37647c6ae99SBarry Smith 
37747c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
37847c6ae99SBarry Smith   dd->Nlocal = x*y*z;
37947c6ae99SBarry Smith   ierr = VecCreateMPIWithArray(comm,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr);
38047c6ae99SBarry Smith   ierr = VecSetBlockSize(global,dof);CHKERRQ(ierr);
38147c6ae99SBarry Smith   dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs);
38247c6ae99SBarry Smith   ierr = VecCreateSeqWithArray(MPI_COMM_SELF,dd->nlocal,0,&local);CHKERRQ(ierr);
38347c6ae99SBarry Smith   ierr = VecSetBlockSize(local,dof);CHKERRQ(ierr);
38447c6ae99SBarry Smith 
38547c6ae99SBarry Smith   /* generate appropriate vector scatters */
38647c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
38747c6ae99SBarry Smith   ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr);
38847c6ae99SBarry Smith   ierr = ISCreateStride(comm,x*y*z,start,1,&to);CHKERRQ(ierr);
38947c6ae99SBarry Smith 
39047c6ae99SBarry Smith   left   = xs - Xs;
39147c6ae99SBarry Smith   bottom = ys - Ys; top = bottom + y;
39247c6ae99SBarry Smith   down   = zs - Zs; up  = down + z;
39347c6ae99SBarry Smith   count  = x*(top-bottom)*(up-down);
39447c6ae99SBarry Smith   ierr = PetscMalloc(count*sizeof(PetscInt)/dof,&idx);CHKERRQ(ierr);
39547c6ae99SBarry Smith   count  = 0;
39647c6ae99SBarry Smith   for (i=down; i<up; i++) {
39747c6ae99SBarry Smith     for (j=bottom; j<top; j++) {
39847c6ae99SBarry Smith       for (k=0; k<x; k += dof) {
39947c6ae99SBarry Smith         idx[count++] = ((left+j*(Xe-Xs))+i*(Xe-Xs)*(Ye-Ys) + k)/dof;
40047c6ae99SBarry Smith       }
40147c6ae99SBarry Smith     }
40247c6ae99SBarry Smith   }
40347c6ae99SBarry Smith 
40447c6ae99SBarry Smith   ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
40547c6ae99SBarry Smith 
40647c6ae99SBarry Smith   ierr = VecScatterCreate(local,from,global,to,&ltog);CHKERRQ(ierr);
40747c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,ltog);CHKERRQ(ierr);
40847c6ae99SBarry Smith   ierr = ISDestroy(from);CHKERRQ(ierr);
40947c6ae99SBarry Smith   ierr = ISDestroy(to);CHKERRQ(ierr);
41047c6ae99SBarry Smith 
41147c6ae99SBarry Smith   /* global to local must include ghost points */
41247c6ae99SBarry Smith   if (stencil_type == DA_STENCIL_BOX) {
41347c6ae99SBarry Smith     ierr = ISCreateStride(comm,(Xe-Xs)*(Ye-Ys)*(Ze-Zs),0,1,&to);CHKERRQ(ierr);
41447c6ae99SBarry Smith   } else {
41547c6ae99SBarry Smith     /* This is way ugly! We need to list the funny cross type region */
41647c6ae99SBarry Smith     /* the bottom chunck */
41747c6ae99SBarry Smith     left   = xs - Xs;
41847c6ae99SBarry Smith     bottom = ys - Ys; top = bottom + y;
41947c6ae99SBarry Smith     down   = zs - Zs;   up  = down + z;
42047c6ae99SBarry Smith     count  = down*(top-bottom)*x + (up-down)*(bottom*x  + (top-bottom)*(Xe-Xs) + (Ye-Ys-top)*x) + (Ze-Zs-up)*(top-bottom)*x;
42147c6ae99SBarry Smith     ierr   = PetscMalloc(count*sizeof(PetscInt)/dof,&idx);CHKERRQ(ierr);
42247c6ae99SBarry Smith     count  = 0;
42347c6ae99SBarry Smith     for (i=0; i<down; i++) {
42447c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
42547c6ae99SBarry Smith         for (k=0; k<x; k += dof) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
42647c6ae99SBarry Smith       }
42747c6ae99SBarry Smith     }
42847c6ae99SBarry Smith     /* the middle piece */
42947c6ae99SBarry Smith     for (i=down; i<up; i++) {
43047c6ae99SBarry Smith       /* front */
43147c6ae99SBarry Smith       for (j=0; j<bottom; j++) {
43247c6ae99SBarry Smith         for (k=0; k<x; k += dof) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
43347c6ae99SBarry Smith       }
43447c6ae99SBarry Smith       /* middle */
43547c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
43647c6ae99SBarry Smith         for (k=0; k<Xe-Xs; k += dof) idx[count++] = j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
43747c6ae99SBarry Smith       }
43847c6ae99SBarry Smith       /* back */
43947c6ae99SBarry Smith       for (j=top; j<Ye-Ys; j++) {
44047c6ae99SBarry Smith         for (k=0; k<x; k += dof) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
44147c6ae99SBarry Smith       }
44247c6ae99SBarry Smith     }
44347c6ae99SBarry Smith     /* the top piece */
44447c6ae99SBarry Smith     for (i=up; i<Ze-Zs; i++) {
44547c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
44647c6ae99SBarry Smith         for (k=0; k<x; k += dof) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
44747c6ae99SBarry Smith       }
44847c6ae99SBarry Smith     }
44947c6ae99SBarry Smith     for (i=0; i<count; i++) idx[i] = idx[i]/dof;
45047c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
45147c6ae99SBarry Smith   }
45247c6ae99SBarry Smith 
45347c6ae99SBarry Smith   /* determine who lies on each side of use stored in    n24 n25 n26
45447c6ae99SBarry Smith                                                          n21 n22 n23
45547c6ae99SBarry Smith                                                          n18 n19 n20
45647c6ae99SBarry Smith 
45747c6ae99SBarry Smith                                                          n15 n16 n17
45847c6ae99SBarry Smith                                                          n12     n14
45947c6ae99SBarry Smith                                                          n9  n10 n11
46047c6ae99SBarry Smith 
46147c6ae99SBarry Smith                                                          n6  n7  n8
46247c6ae99SBarry Smith                                                          n3  n4  n5
46347c6ae99SBarry Smith                                                          n0  n1  n2
46447c6ae99SBarry Smith   */
46547c6ae99SBarry Smith 
46647c6ae99SBarry Smith   /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
46747c6ae99SBarry Smith 
46847c6ae99SBarry Smith   /* Assume Nodes are Internal to the Cube */
46947c6ae99SBarry Smith 
47047c6ae99SBarry Smith   n0  = rank - m*n - m - 1;
47147c6ae99SBarry Smith   n1  = rank - m*n - m;
47247c6ae99SBarry Smith   n2  = rank - m*n - m + 1;
47347c6ae99SBarry Smith   n3  = rank - m*n -1;
47447c6ae99SBarry Smith   n4  = rank - m*n;
47547c6ae99SBarry Smith   n5  = rank - m*n + 1;
47647c6ae99SBarry Smith   n6  = rank - m*n + m - 1;
47747c6ae99SBarry Smith   n7  = rank - m*n + m;
47847c6ae99SBarry Smith   n8  = rank - m*n + m + 1;
47947c6ae99SBarry Smith 
48047c6ae99SBarry Smith   n9  = rank - m - 1;
48147c6ae99SBarry Smith   n10 = rank - m;
48247c6ae99SBarry Smith   n11 = rank - m + 1;
48347c6ae99SBarry Smith   n12 = rank - 1;
48447c6ae99SBarry Smith   n14 = rank + 1;
48547c6ae99SBarry Smith   n15 = rank + m - 1;
48647c6ae99SBarry Smith   n16 = rank + m;
48747c6ae99SBarry Smith   n17 = rank + m + 1;
48847c6ae99SBarry Smith 
48947c6ae99SBarry Smith   n18 = rank + m*n - m - 1;
49047c6ae99SBarry Smith   n19 = rank + m*n - m;
49147c6ae99SBarry Smith   n20 = rank + m*n - m + 1;
49247c6ae99SBarry Smith   n21 = rank + m*n - 1;
49347c6ae99SBarry Smith   n22 = rank + m*n;
49447c6ae99SBarry Smith   n23 = rank + m*n + 1;
49547c6ae99SBarry Smith   n24 = rank + m*n + m - 1;
49647c6ae99SBarry Smith   n25 = rank + m*n + m;
49747c6ae99SBarry Smith   n26 = rank + m*n + m + 1;
49847c6ae99SBarry Smith 
49947c6ae99SBarry Smith   /* Assume Pieces are on Faces of Cube */
50047c6ae99SBarry Smith 
50147c6ae99SBarry Smith   if (xs == 0) { /* First assume not corner or edge */
50247c6ae99SBarry Smith     n0  = rank       -1 - (m*n);
50347c6ae99SBarry Smith     n3  = rank + m   -1 - (m*n);
50447c6ae99SBarry Smith     n6  = rank + 2*m -1 - (m*n);
50547c6ae99SBarry Smith     n9  = rank       -1;
50647c6ae99SBarry Smith     n12 = rank + m   -1;
50747c6ae99SBarry Smith     n15 = rank + 2*m -1;
50847c6ae99SBarry Smith     n18 = rank       -1 + (m*n);
50947c6ae99SBarry Smith     n21 = rank + m   -1 + (m*n);
51047c6ae99SBarry Smith     n24 = rank + 2*m -1 + (m*n);
51147c6ae99SBarry Smith    }
51247c6ae99SBarry Smith 
51347c6ae99SBarry Smith   if (xe == M*dof) { /* First assume not corner or edge */
51447c6ae99SBarry Smith     n2  = rank -2*m +1 - (m*n);
51547c6ae99SBarry Smith     n5  = rank - m  +1 - (m*n);
51647c6ae99SBarry Smith     n8  = rank      +1 - (m*n);
51747c6ae99SBarry Smith     n11 = rank -2*m +1;
51847c6ae99SBarry Smith     n14 = rank - m  +1;
51947c6ae99SBarry Smith     n17 = rank      +1;
52047c6ae99SBarry Smith     n20 = rank -2*m +1 + (m*n);
52147c6ae99SBarry Smith     n23 = rank - m  +1 + (m*n);
52247c6ae99SBarry Smith     n26 = rank      +1 + (m*n);
52347c6ae99SBarry Smith   }
52447c6ae99SBarry Smith 
52547c6ae99SBarry Smith   if (ys==0) { /* First assume not corner or edge */
52647c6ae99SBarry Smith     n0  = rank + m * (n-1) -1 - (m*n);
52747c6ae99SBarry Smith     n1  = rank + m * (n-1)    - (m*n);
52847c6ae99SBarry Smith     n2  = rank + m * (n-1) +1 - (m*n);
52947c6ae99SBarry Smith     n9  = rank + m * (n-1) -1;
53047c6ae99SBarry Smith     n10 = rank + m * (n-1);
53147c6ae99SBarry Smith     n11 = rank + m * (n-1) +1;
53247c6ae99SBarry Smith     n18 = rank + m * (n-1) -1 + (m*n);
53347c6ae99SBarry Smith     n19 = rank + m * (n-1)    + (m*n);
53447c6ae99SBarry Smith     n20 = rank + m * (n-1) +1 + (m*n);
53547c6ae99SBarry Smith   }
53647c6ae99SBarry Smith 
53747c6ae99SBarry Smith   if (ye == N) { /* First assume not corner or edge */
53847c6ae99SBarry Smith     n6  = rank - m * (n-1) -1 - (m*n);
53947c6ae99SBarry Smith     n7  = rank - m * (n-1)    - (m*n);
54047c6ae99SBarry Smith     n8  = rank - m * (n-1) +1 - (m*n);
54147c6ae99SBarry Smith     n15 = rank - m * (n-1) -1;
54247c6ae99SBarry Smith     n16 = rank - m * (n-1);
54347c6ae99SBarry Smith     n17 = rank - m * (n-1) +1;
54447c6ae99SBarry Smith     n24 = rank - m * (n-1) -1 + (m*n);
54547c6ae99SBarry Smith     n25 = rank - m * (n-1)    + (m*n);
54647c6ae99SBarry Smith     n26 = rank - m * (n-1) +1 + (m*n);
54747c6ae99SBarry Smith   }
54847c6ae99SBarry Smith 
54947c6ae99SBarry Smith   if (zs == 0) { /* First assume not corner or edge */
55047c6ae99SBarry Smith     n0 = size - (m*n) + rank - m - 1;
55147c6ae99SBarry Smith     n1 = size - (m*n) + rank - m;
55247c6ae99SBarry Smith     n2 = size - (m*n) + rank - m + 1;
55347c6ae99SBarry Smith     n3 = size - (m*n) + rank - 1;
55447c6ae99SBarry Smith     n4 = size - (m*n) + rank;
55547c6ae99SBarry Smith     n5 = size - (m*n) + rank + 1;
55647c6ae99SBarry Smith     n6 = size - (m*n) + rank + m - 1;
55747c6ae99SBarry Smith     n7 = size - (m*n) + rank + m ;
55847c6ae99SBarry Smith     n8 = size - (m*n) + rank + m + 1;
55947c6ae99SBarry Smith   }
56047c6ae99SBarry Smith 
56147c6ae99SBarry Smith   if (ze == P) { /* First assume not corner or edge */
56247c6ae99SBarry Smith     n18 = (m*n) - (size-rank) - m - 1;
56347c6ae99SBarry Smith     n19 = (m*n) - (size-rank) - m;
56447c6ae99SBarry Smith     n20 = (m*n) - (size-rank) - m + 1;
56547c6ae99SBarry Smith     n21 = (m*n) - (size-rank) - 1;
56647c6ae99SBarry Smith     n22 = (m*n) - (size-rank);
56747c6ae99SBarry Smith     n23 = (m*n) - (size-rank) + 1;
56847c6ae99SBarry Smith     n24 = (m*n) - (size-rank) + m - 1;
56947c6ae99SBarry Smith     n25 = (m*n) - (size-rank) + m;
57047c6ae99SBarry Smith     n26 = (m*n) - (size-rank) + m + 1;
57147c6ae99SBarry Smith   }
57247c6ae99SBarry Smith 
57347c6ae99SBarry Smith   if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */
57447c6ae99SBarry Smith     n0 = size - m*n + rank + m-1 - m;
57547c6ae99SBarry Smith     n3 = size - m*n + rank + m-1;
57647c6ae99SBarry Smith     n6 = size - m*n + rank + m-1 + m;
57747c6ae99SBarry Smith   }
57847c6ae99SBarry Smith 
57947c6ae99SBarry Smith   if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */
58047c6ae99SBarry Smith     n18 = m*n - (size - rank) + m-1 - m;
58147c6ae99SBarry Smith     n21 = m*n - (size - rank) + m-1;
58247c6ae99SBarry Smith     n24 = m*n - (size - rank) + m-1 + m;
58347c6ae99SBarry Smith   }
58447c6ae99SBarry Smith 
58547c6ae99SBarry Smith   if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */
58647c6ae99SBarry Smith     n0  = rank + m*n -1 - m*n;
58747c6ae99SBarry Smith     n9  = rank + m*n -1;
58847c6ae99SBarry Smith     n18 = rank + m*n -1 + m*n;
58947c6ae99SBarry Smith   }
59047c6ae99SBarry Smith 
59147c6ae99SBarry Smith   if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */
59247c6ae99SBarry Smith     n6  = rank - m*(n-1) + m-1 - m*n;
59347c6ae99SBarry Smith     n15 = rank - m*(n-1) + m-1;
59447c6ae99SBarry Smith     n24 = rank - m*(n-1) + m-1 + m*n;
59547c6ae99SBarry Smith   }
59647c6ae99SBarry Smith 
59747c6ae99SBarry Smith   if ((xe==M*dof) && (zs==0)) { /* Assume an edge, not corner */
59847c6ae99SBarry Smith     n2 = size - (m*n-rank) - (m-1) - m;
59947c6ae99SBarry Smith     n5 = size - (m*n-rank) - (m-1);
60047c6ae99SBarry Smith     n8 = size - (m*n-rank) - (m-1) + m;
60147c6ae99SBarry Smith   }
60247c6ae99SBarry Smith 
60347c6ae99SBarry Smith   if ((xe==M*dof) && (ze==P)) { /* Assume an edge, not corner */
60447c6ae99SBarry Smith     n20 = m*n - (size - rank) - (m-1) - m;
60547c6ae99SBarry Smith     n23 = m*n - (size - rank) - (m-1);
60647c6ae99SBarry Smith     n26 = m*n - (size - rank) - (m-1) + m;
60747c6ae99SBarry Smith   }
60847c6ae99SBarry Smith 
60947c6ae99SBarry Smith   if ((xe==M*dof) && (ys==0)) { /* Assume an edge, not corner */
61047c6ae99SBarry Smith     n2  = rank + m*(n-1) - (m-1) - m*n;
61147c6ae99SBarry Smith     n11 = rank + m*(n-1) - (m-1);
61247c6ae99SBarry Smith     n20 = rank + m*(n-1) - (m-1) + m*n;
61347c6ae99SBarry Smith   }
61447c6ae99SBarry Smith 
61547c6ae99SBarry Smith   if ((xe==M*dof) && (ye==N)) { /* Assume an edge, not corner */
61647c6ae99SBarry Smith     n8  = rank - m*n +1 - m*n;
61747c6ae99SBarry Smith     n17 = rank - m*n +1;
61847c6ae99SBarry Smith     n26 = rank - m*n +1 + m*n;
61947c6ae99SBarry Smith   }
62047c6ae99SBarry Smith 
62147c6ae99SBarry Smith   if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */
62247c6ae99SBarry Smith     n0 = size - m + rank -1;
62347c6ae99SBarry Smith     n1 = size - m + rank;
62447c6ae99SBarry Smith     n2 = size - m + rank +1;
62547c6ae99SBarry Smith   }
62647c6ae99SBarry Smith 
62747c6ae99SBarry Smith   if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */
62847c6ae99SBarry Smith     n18 = m*n - (size - rank) + m*(n-1) -1;
62947c6ae99SBarry Smith     n19 = m*n - (size - rank) + m*(n-1);
63047c6ae99SBarry Smith     n20 = m*n - (size - rank) + m*(n-1) +1;
63147c6ae99SBarry Smith   }
63247c6ae99SBarry Smith 
63347c6ae99SBarry Smith   if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */
63447c6ae99SBarry Smith     n6 = size - (m*n-rank) - m * (n-1) -1;
63547c6ae99SBarry Smith     n7 = size - (m*n-rank) - m * (n-1);
63647c6ae99SBarry Smith     n8 = size - (m*n-rank) - m * (n-1) +1;
63747c6ae99SBarry Smith   }
63847c6ae99SBarry Smith 
63947c6ae99SBarry Smith   if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */
64047c6ae99SBarry Smith     n24 = rank - (size-m) -1;
64147c6ae99SBarry Smith     n25 = rank - (size-m);
64247c6ae99SBarry Smith     n26 = rank - (size-m) +1;
64347c6ae99SBarry Smith   }
64447c6ae99SBarry Smith 
64547c6ae99SBarry Smith   /* Check for Corners */
64647c6ae99SBarry Smith   if ((xs==0)   && (ys==0) && (zs==0)) { n0  = size -1;}
64747c6ae99SBarry Smith   if ((xs==0)   && (ys==0) && (ze==P)) { n18 = m*n-1;}
64847c6ae99SBarry Smith   if ((xs==0)   && (ye==N) && (zs==0)) { n6  = (size-1)-m*(n-1);}
64947c6ae99SBarry Smith   if ((xs==0)   && (ye==N) && (ze==P)) { n24 = m-1;}
65047c6ae99SBarry Smith   if ((xe==M*dof) && (ys==0) && (zs==0)) { n2  = size-m;}
65147c6ae99SBarry Smith   if ((xe==M*dof) && (ys==0) && (ze==P)) { n20 = m*n-m;}
65247c6ae99SBarry Smith   if ((xe==M*dof) && (ye==N) && (zs==0)) { n8  = size-m*n;}
65347c6ae99SBarry Smith   if ((xe==M*dof) && (ye==N) && (ze==P)) { n26 = 0;}
65447c6ae99SBarry Smith 
65547c6ae99SBarry Smith   /* Check for when not X,Y, and Z Periodic */
65647c6ae99SBarry Smith 
65747c6ae99SBarry Smith   /* If not X periodic */
65847c6ae99SBarry Smith   if ((wrap != DA_XPERIODIC)  && (wrap != DA_XYPERIODIC) &&
65947c6ae99SBarry Smith      (wrap != DA_XZPERIODIC) && (wrap != DA_XYZPERIODIC)) {
66047c6ae99SBarry Smith     if (xs==0)   {n0  = n3  = n6  = n9  = n12 = n15 = n18 = n21 = n24 = -2;}
66147c6ae99SBarry Smith     if (xe==M*dof) {n2  = n5  = n8  = n11 = n14 = n17 = n20 = n23 = n26 = -2;}
66247c6ae99SBarry Smith   }
66347c6ae99SBarry Smith 
66447c6ae99SBarry Smith   /* If not Y periodic */
66547c6ae99SBarry Smith   if ((wrap != DA_YPERIODIC)  && (wrap != DA_XYPERIODIC) &&
66647c6ae99SBarry Smith       (wrap != DA_YZPERIODIC) && (wrap != DA_XYZPERIODIC)) {
66747c6ae99SBarry Smith     if (ys==0)   {n0  = n1  = n2  = n9  = n10 = n11 = n18 = n19 = n20 = -2;}
66847c6ae99SBarry Smith     if (ye==N)   {n6  = n7  = n8  = n15 = n16 = n17 = n24 = n25 = n26 = -2;}
66947c6ae99SBarry Smith   }
67047c6ae99SBarry Smith 
67147c6ae99SBarry Smith   /* If not Z periodic */
67247c6ae99SBarry Smith   if ((wrap != DA_ZPERIODIC)  && (wrap != DA_XZPERIODIC) &&
67347c6ae99SBarry Smith       (wrap != DA_YZPERIODIC) && (wrap != DA_XYZPERIODIC)) {
67447c6ae99SBarry Smith     if (zs==0)   {n0  = n1  = n2  = n3  = n4  = n5  = n6  = n7  = n8  = -2;}
67547c6ae99SBarry Smith     if (ze==P)   {n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;}
67647c6ae99SBarry Smith   }
67747c6ae99SBarry Smith 
67847c6ae99SBarry Smith   ierr = PetscMalloc(27*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr);
67947c6ae99SBarry Smith   dd->neighbors[0] = n0;
68047c6ae99SBarry Smith   dd->neighbors[1] = n1;
68147c6ae99SBarry Smith   dd->neighbors[2] = n2;
68247c6ae99SBarry Smith   dd->neighbors[3] = n3;
68347c6ae99SBarry Smith   dd->neighbors[4] = n4;
68447c6ae99SBarry Smith   dd->neighbors[5] = n5;
68547c6ae99SBarry Smith   dd->neighbors[6] = n6;
68647c6ae99SBarry Smith   dd->neighbors[7] = n7;
68747c6ae99SBarry Smith   dd->neighbors[8] = n8;
68847c6ae99SBarry Smith   dd->neighbors[9] = n9;
68947c6ae99SBarry Smith   dd->neighbors[10] = n10;
69047c6ae99SBarry Smith   dd->neighbors[11] = n11;
69147c6ae99SBarry Smith   dd->neighbors[12] = n12;
69247c6ae99SBarry Smith   dd->neighbors[13] = rank;
69347c6ae99SBarry Smith   dd->neighbors[14] = n14;
69447c6ae99SBarry Smith   dd->neighbors[15] = n15;
69547c6ae99SBarry Smith   dd->neighbors[16] = n16;
69647c6ae99SBarry Smith   dd->neighbors[17] = n17;
69747c6ae99SBarry Smith   dd->neighbors[18] = n18;
69847c6ae99SBarry Smith   dd->neighbors[19] = n19;
69947c6ae99SBarry Smith   dd->neighbors[20] = n20;
70047c6ae99SBarry Smith   dd->neighbors[21] = n21;
70147c6ae99SBarry Smith   dd->neighbors[22] = n22;
70247c6ae99SBarry Smith   dd->neighbors[23] = n23;
70347c6ae99SBarry Smith   dd->neighbors[24] = n24;
70447c6ae99SBarry Smith   dd->neighbors[25] = n25;
70547c6ae99SBarry Smith   dd->neighbors[26] = n26;
70647c6ae99SBarry Smith 
70747c6ae99SBarry Smith   /* If star stencil then delete the corner neighbors */
70847c6ae99SBarry Smith   if (stencil_type == DA_STENCIL_STAR) {
70947c6ae99SBarry Smith      /* save information about corner neighbors */
71047c6ae99SBarry Smith      sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7;
71147c6ae99SBarry Smith      sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18;
71247c6ae99SBarry Smith      sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25;
71347c6ae99SBarry Smith      sn26 = n26;
71447c6ae99SBarry Smith      n0  = n1  = n2  = n3  = n5  = n6  = n7  = n8  = n9  = n11 =
71547c6ae99SBarry Smith      n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1;
71647c6ae99SBarry Smith   }
71747c6ae99SBarry Smith 
71847c6ae99SBarry Smith 
71947c6ae99SBarry Smith   ierr = PetscMalloc((Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
72047c6ae99SBarry Smith   ierr = PetscLogObjectMemory(da,(Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt));CHKERRQ(ierr);
72147c6ae99SBarry Smith 
72247c6ae99SBarry Smith   nn = 0;
72347c6ae99SBarry Smith 
72447c6ae99SBarry Smith   /* Bottom Level */
72547c6ae99SBarry Smith   for (k=0; k<s_z; k++) {
72647c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
72747c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
72847c6ae99SBarry Smith         x_t = lx[n0 % m]*dof;
72947c6ae99SBarry Smith         y_t = ly[(n0 % (m*n))/m];
73047c6ae99SBarry Smith         z_t = lz[n0 / (m*n)];
73147c6ae99SBarry Smith         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;
73247c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
73347c6ae99SBarry Smith       }
73447c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
73547c6ae99SBarry Smith         x_t = x;
73647c6ae99SBarry Smith         y_t = ly[(n1 % (m*n))/m];
73747c6ae99SBarry Smith         z_t = lz[n1 / (m*n)];
73847c6ae99SBarry Smith         s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
73947c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
74047c6ae99SBarry Smith       }
74147c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
74247c6ae99SBarry Smith         x_t = lx[n2 % m]*dof;
74347c6ae99SBarry Smith         y_t = ly[(n2 % (m*n))/m];
74447c6ae99SBarry Smith         z_t = lz[n2 / (m*n)];
74547c6ae99SBarry Smith         s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
74647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
74747c6ae99SBarry Smith       }
74847c6ae99SBarry Smith     }
74947c6ae99SBarry Smith 
75047c6ae99SBarry Smith     for (i=0; i<y; i++) {
75147c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
75247c6ae99SBarry Smith         x_t = lx[n3 % m]*dof;
75347c6ae99SBarry Smith         y_t = y;
75447c6ae99SBarry Smith         z_t = lz[n3 / (m*n)];
75547c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
75647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
75747c6ae99SBarry Smith       }
75847c6ae99SBarry Smith 
75947c6ae99SBarry Smith       if (n4 >= 0) { /* middle */
76047c6ae99SBarry Smith         x_t = x;
76147c6ae99SBarry Smith         y_t = y;
76247c6ae99SBarry Smith         z_t = lz[n4 / (m*n)];
76347c6ae99SBarry Smith         s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
76447c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
76547c6ae99SBarry Smith       }
76647c6ae99SBarry Smith 
76747c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
76847c6ae99SBarry Smith         x_t = lx[n5 % m]*dof;
76947c6ae99SBarry Smith         y_t = y;
77047c6ae99SBarry Smith         z_t = lz[n5 / (m*n)];
77147c6ae99SBarry Smith         s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
77247c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
77347c6ae99SBarry Smith       }
77447c6ae99SBarry Smith     }
77547c6ae99SBarry Smith 
77647c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
77747c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
77847c6ae99SBarry Smith         x_t = lx[n6 % m]*dof;
77947c6ae99SBarry Smith         y_t = ly[(n6 % (m*n))/m];
78047c6ae99SBarry Smith         z_t = lz[n6 / (m*n)];
78147c6ae99SBarry Smith         s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
78247c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
78347c6ae99SBarry Smith       }
78447c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
78547c6ae99SBarry Smith         x_t = x;
78647c6ae99SBarry Smith         y_t = ly[(n7 % (m*n))/m];
78747c6ae99SBarry Smith         z_t = lz[n7 / (m*n)];
78847c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
78947c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
79047c6ae99SBarry Smith       }
79147c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
79247c6ae99SBarry Smith         x_t = lx[n8 % m]*dof;
79347c6ae99SBarry Smith         y_t = ly[(n8 % (m*n))/m];
79447c6ae99SBarry Smith         z_t = lz[n8 / (m*n)];
79547c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
79647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
79747c6ae99SBarry Smith       }
79847c6ae99SBarry Smith     }
79947c6ae99SBarry Smith   }
80047c6ae99SBarry Smith 
80147c6ae99SBarry Smith   /* Middle Level */
80247c6ae99SBarry Smith   for (k=0; k<z; k++) {
80347c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
80447c6ae99SBarry Smith       if (n9 >= 0) { /* left below */
80547c6ae99SBarry Smith         x_t = lx[n9 % m]*dof;
80647c6ae99SBarry Smith         y_t = ly[(n9 % (m*n))/m];
80747c6ae99SBarry Smith         /* z_t = z; */
80847c6ae99SBarry Smith         s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
80947c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
81047c6ae99SBarry Smith       }
81147c6ae99SBarry Smith       if (n10 >= 0) { /* directly below */
81247c6ae99SBarry Smith         x_t = x;
81347c6ae99SBarry Smith         y_t = ly[(n10 % (m*n))/m];
81447c6ae99SBarry Smith         /* z_t = z; */
81547c6ae99SBarry Smith         s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
81647c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
81747c6ae99SBarry Smith       }
81847c6ae99SBarry Smith       if (n11 >= 0) { /* right below */
81947c6ae99SBarry Smith         x_t = lx[n11 % m]*dof;
82047c6ae99SBarry Smith         y_t = ly[(n11 % (m*n))/m];
82147c6ae99SBarry Smith         /* z_t = z; */
82247c6ae99SBarry Smith         s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
82347c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
82447c6ae99SBarry Smith       }
82547c6ae99SBarry Smith     }
82647c6ae99SBarry Smith 
82747c6ae99SBarry Smith     for (i=0; i<y; i++) {
82847c6ae99SBarry Smith       if (n12 >= 0) { /* directly left */
82947c6ae99SBarry Smith         x_t = lx[n12 % m]*dof;
83047c6ae99SBarry Smith         y_t = y;
83147c6ae99SBarry Smith         /* z_t = z; */
83247c6ae99SBarry Smith         s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
83347c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
83447c6ae99SBarry Smith       }
83547c6ae99SBarry Smith 
83647c6ae99SBarry Smith       /* Interior */
83747c6ae99SBarry Smith       s_t = bases[rank] + i*x + k*x*y;
83847c6ae99SBarry Smith       for (j=0; j<x; j++) { idx[nn++] = s_t++;}
83947c6ae99SBarry Smith 
84047c6ae99SBarry Smith       if (n14 >= 0) { /* directly right */
84147c6ae99SBarry Smith         x_t = lx[n14 % m]*dof;
84247c6ae99SBarry Smith         y_t = y;
84347c6ae99SBarry Smith         /* z_t = z; */
84447c6ae99SBarry Smith         s_t = bases[n14] + i*x_t + k*x_t*y_t;
84547c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
84647c6ae99SBarry Smith       }
84747c6ae99SBarry Smith     }
84847c6ae99SBarry Smith 
84947c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
85047c6ae99SBarry Smith       if (n15 >= 0) { /* left above */
85147c6ae99SBarry Smith         x_t = lx[n15 % m]*dof;
85247c6ae99SBarry Smith         y_t = ly[(n15 % (m*n))/m];
85347c6ae99SBarry Smith         /* z_t = z; */
85447c6ae99SBarry Smith         s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
85547c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
85647c6ae99SBarry Smith       }
85747c6ae99SBarry Smith       if (n16 >= 0) { /* directly above */
85847c6ae99SBarry Smith         x_t = x;
85947c6ae99SBarry Smith         y_t = ly[(n16 % (m*n))/m];
86047c6ae99SBarry Smith         /* z_t = z; */
86147c6ae99SBarry Smith         s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
86247c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
86347c6ae99SBarry Smith       }
86447c6ae99SBarry Smith       if (n17 >= 0) { /* right above */
86547c6ae99SBarry Smith         x_t = lx[n17 % m]*dof;
86647c6ae99SBarry Smith         y_t = ly[(n17 % (m*n))/m];
86747c6ae99SBarry Smith         /* z_t = z; */
86847c6ae99SBarry Smith         s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
86947c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
87047c6ae99SBarry Smith       }
87147c6ae99SBarry Smith     }
87247c6ae99SBarry Smith   }
87347c6ae99SBarry Smith 
87447c6ae99SBarry Smith   /* Upper Level */
87547c6ae99SBarry Smith   for (k=0; k<s_z; k++) {
87647c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
87747c6ae99SBarry Smith       if (n18 >= 0) { /* left below */
87847c6ae99SBarry Smith         x_t = lx[n18 % m]*dof;
87947c6ae99SBarry Smith         y_t = ly[(n18 % (m*n))/m];
88047c6ae99SBarry Smith         /* z_t = lz[n18 / (m*n)]; */
88147c6ae99SBarry Smith         s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
88247c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
88347c6ae99SBarry Smith       }
88447c6ae99SBarry Smith       if (n19 >= 0) { /* directly below */
88547c6ae99SBarry Smith         x_t = x;
88647c6ae99SBarry Smith         y_t = ly[(n19 % (m*n))/m];
88747c6ae99SBarry Smith         /* z_t = lz[n19 / (m*n)]; */
88847c6ae99SBarry Smith         s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
88947c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
89047c6ae99SBarry Smith       }
89147c6ae99SBarry Smith       if (n20 >= 0) { /* right below */
89247c6ae99SBarry Smith         x_t = lx[n20 % m]*dof;
89347c6ae99SBarry Smith         y_t = ly[(n20 % (m*n))/m];
89447c6ae99SBarry Smith         /* z_t = lz[n20 / (m*n)]; */
89547c6ae99SBarry Smith         s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
89647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
89747c6ae99SBarry Smith       }
89847c6ae99SBarry Smith     }
89947c6ae99SBarry Smith 
90047c6ae99SBarry Smith     for (i=0; i<y; i++) {
90147c6ae99SBarry Smith       if (n21 >= 0) { /* directly left */
90247c6ae99SBarry Smith         x_t = lx[n21 % m]*dof;
90347c6ae99SBarry Smith         y_t = y;
90447c6ae99SBarry Smith         /* z_t = lz[n21 / (m*n)]; */
90547c6ae99SBarry Smith         s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
90647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
90747c6ae99SBarry Smith       }
90847c6ae99SBarry Smith 
90947c6ae99SBarry Smith       if (n22 >= 0) { /* middle */
91047c6ae99SBarry Smith         x_t = x;
91147c6ae99SBarry Smith         y_t = y;
91247c6ae99SBarry Smith         /* z_t = lz[n22 / (m*n)]; */
91347c6ae99SBarry Smith         s_t = bases[n22] + i*x_t + k*x_t*y_t;
91447c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
91547c6ae99SBarry Smith       }
91647c6ae99SBarry Smith 
91747c6ae99SBarry Smith       if (n23 >= 0) { /* directly right */
91847c6ae99SBarry Smith         x_t = lx[n23 % m]*dof;
91947c6ae99SBarry Smith         y_t = y;
92047c6ae99SBarry Smith         /* z_t = lz[n23 / (m*n)]; */
92147c6ae99SBarry Smith         s_t = bases[n23] + i*x_t + k*x_t*y_t;
92247c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
92347c6ae99SBarry Smith       }
92447c6ae99SBarry Smith     }
92547c6ae99SBarry Smith 
92647c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
92747c6ae99SBarry Smith       if (n24 >= 0) { /* left above */
92847c6ae99SBarry Smith         x_t = lx[n24 % m]*dof;
92947c6ae99SBarry Smith         y_t = ly[(n24 % (m*n))/m];
93047c6ae99SBarry Smith         /* z_t = lz[n24 / (m*n)]; */
93147c6ae99SBarry Smith         s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
93247c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
93347c6ae99SBarry Smith       }
93447c6ae99SBarry Smith       if (n25 >= 0) { /* directly above */
93547c6ae99SBarry Smith         x_t = x;
93647c6ae99SBarry Smith         y_t = ly[(n25 % (m*n))/m];
93747c6ae99SBarry Smith         /* z_t = lz[n25 / (m*n)]; */
93847c6ae99SBarry Smith         s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
93947c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
94047c6ae99SBarry Smith       }
94147c6ae99SBarry Smith       if (n26 >= 0) { /* right above */
94247c6ae99SBarry Smith         x_t = lx[n26 % m]*dof;
94347c6ae99SBarry Smith         y_t = ly[(n26 % (m*n))/m];
94447c6ae99SBarry Smith         /* z_t = lz[n26 / (m*n)]; */
94547c6ae99SBarry Smith         s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
94647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
94747c6ae99SBarry Smith       }
94847c6ae99SBarry Smith     }
94947c6ae99SBarry Smith   }
95047c6ae99SBarry Smith   base = bases[rank];
95147c6ae99SBarry Smith   {
95247c6ae99SBarry Smith     PetscInt nnn = nn/dof,*iidx;
95347c6ae99SBarry Smith     ierr = PetscMalloc(nnn*sizeof(PetscInt),&iidx);CHKERRQ(ierr);
95447c6ae99SBarry Smith     for (i=0; i<nnn; i++) {
95547c6ae99SBarry Smith       iidx[i] = idx[dof*i]/dof;
95647c6ae99SBarry Smith     }
95747c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,nnn,iidx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
95847c6ae99SBarry Smith   }
95947c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
96047c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr);
96147c6ae99SBarry Smith   ierr = ISDestroy(from);CHKERRQ(ierr);
96247c6ae99SBarry Smith   ierr = ISDestroy(to);CHKERRQ(ierr);
96347c6ae99SBarry Smith 
96447c6ae99SBarry Smith   dd->m  = m;  dd->n  = n;  dd->p  = p;
96547c6ae99SBarry Smith   dd->xs = xs; dd->xe = xe; dd->ys = ys; dd->ye = ye; dd->zs = zs; dd->ze = ze;
96647c6ae99SBarry Smith   dd->Xs = Xs; dd->Xe = Xe; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = Zs; dd->Ze = Ze;
96747c6ae99SBarry Smith 
96847c6ae99SBarry Smith   ierr = VecDestroy(local);CHKERRQ(ierr);
96947c6ae99SBarry Smith   ierr = VecDestroy(global);CHKERRQ(ierr);
97047c6ae99SBarry Smith 
97147c6ae99SBarry Smith   if (stencil_type == DA_STENCIL_STAR) {
97247c6ae99SBarry Smith     /*
97347c6ae99SBarry Smith         Recompute the local to global mappings, this time keeping the
97447c6ae99SBarry Smith       information about the cross corner processor numbers.
97547c6ae99SBarry Smith     */
97647c6ae99SBarry Smith     n0  = sn0;  n1  = sn1;  n2  = sn2;  n3  = sn3;  n5  = sn5;  n6  = sn6; n7 = sn7;
97747c6ae99SBarry Smith     n8  = sn8;  n9  = sn9;  n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18;
97847c6ae99SBarry Smith     n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25;
97947c6ae99SBarry Smith     n26 = sn26;
98047c6ae99SBarry Smith 
98147c6ae99SBarry Smith     nn = 0;
98247c6ae99SBarry Smith 
98347c6ae99SBarry Smith     /* Bottom Level */
98447c6ae99SBarry Smith     for (k=0; k<s_z; k++) {
98547c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
98647c6ae99SBarry Smith         if (n0 >= 0) { /* left below */
98747c6ae99SBarry Smith           x_t = lx[n0 % m]*dof;
98847c6ae99SBarry Smith           y_t = ly[(n0 % (m*n))/m];
98947c6ae99SBarry Smith           z_t = lz[n0 / (m*n)];
99047c6ae99SBarry Smith           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;
99147c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
99247c6ae99SBarry Smith         }
99347c6ae99SBarry Smith         if (n1 >= 0) { /* directly below */
99447c6ae99SBarry Smith           x_t = x;
99547c6ae99SBarry Smith           y_t = ly[(n1 % (m*n))/m];
99647c6ae99SBarry Smith           z_t = lz[n1 / (m*n)];
99747c6ae99SBarry Smith           s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
99847c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
99947c6ae99SBarry Smith         }
100047c6ae99SBarry Smith         if (n2 >= 0) { /* right below */
100147c6ae99SBarry Smith           x_t = lx[n2 % m]*dof;
100247c6ae99SBarry Smith           y_t = ly[(n2 % (m*n))/m];
100347c6ae99SBarry Smith           z_t = lz[n2 / (m*n)];
100447c6ae99SBarry Smith           s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
100547c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
100647c6ae99SBarry Smith         }
100747c6ae99SBarry Smith       }
100847c6ae99SBarry Smith 
100947c6ae99SBarry Smith       for (i=0; i<y; i++) {
101047c6ae99SBarry Smith         if (n3 >= 0) { /* directly left */
101147c6ae99SBarry Smith           x_t = lx[n3 % m]*dof;
101247c6ae99SBarry Smith           y_t = y;
101347c6ae99SBarry Smith           z_t = lz[n3 / (m*n)];
101447c6ae99SBarry Smith           s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
101547c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
101647c6ae99SBarry Smith         }
101747c6ae99SBarry Smith 
101847c6ae99SBarry Smith         if (n4 >= 0) { /* middle */
101947c6ae99SBarry Smith           x_t = x;
102047c6ae99SBarry Smith           y_t = y;
102147c6ae99SBarry Smith           z_t = lz[n4 / (m*n)];
102247c6ae99SBarry Smith           s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
102347c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
102447c6ae99SBarry Smith         }
102547c6ae99SBarry Smith 
102647c6ae99SBarry Smith         if (n5 >= 0) { /* directly right */
102747c6ae99SBarry Smith           x_t = lx[n5 % m]*dof;
102847c6ae99SBarry Smith           y_t = y;
102947c6ae99SBarry Smith           z_t = lz[n5 / (m*n)];
103047c6ae99SBarry Smith           s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
103147c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
103247c6ae99SBarry Smith         }
103347c6ae99SBarry Smith       }
103447c6ae99SBarry Smith 
103547c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
103647c6ae99SBarry Smith         if (n6 >= 0) { /* left above */
103747c6ae99SBarry Smith           x_t = lx[n6 % m]*dof;
103847c6ae99SBarry Smith           y_t = ly[(n6 % (m*n))/m];
103947c6ae99SBarry Smith           z_t = lz[n6 / (m*n)];
104047c6ae99SBarry Smith           s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
104147c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
104247c6ae99SBarry Smith         }
104347c6ae99SBarry Smith         if (n7 >= 0) { /* directly above */
104447c6ae99SBarry Smith           x_t = x;
104547c6ae99SBarry Smith           y_t = ly[(n7 % (m*n))/m];
104647c6ae99SBarry Smith           z_t = lz[n7 / (m*n)];
104747c6ae99SBarry Smith           s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
104847c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
104947c6ae99SBarry Smith         }
105047c6ae99SBarry Smith         if (n8 >= 0) { /* right above */
105147c6ae99SBarry Smith           x_t = lx[n8 % m]*dof;
105247c6ae99SBarry Smith           y_t = ly[(n8 % (m*n))/m];
105347c6ae99SBarry Smith           z_t = lz[n8 / (m*n)];
105447c6ae99SBarry Smith           s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
105547c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
105647c6ae99SBarry Smith         }
105747c6ae99SBarry Smith       }
105847c6ae99SBarry Smith     }
105947c6ae99SBarry Smith 
106047c6ae99SBarry Smith     /* Middle Level */
106147c6ae99SBarry Smith     for (k=0; k<z; k++) {
106247c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
106347c6ae99SBarry Smith         if (n9 >= 0) { /* left below */
106447c6ae99SBarry Smith           x_t = lx[n9 % m]*dof;
106547c6ae99SBarry Smith           y_t = ly[(n9 % (m*n))/m];
106647c6ae99SBarry Smith           /* z_t = z; */
106747c6ae99SBarry Smith           s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
106847c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
106947c6ae99SBarry Smith         }
107047c6ae99SBarry Smith         if (n10 >= 0) { /* directly below */
107147c6ae99SBarry Smith           x_t = x;
107247c6ae99SBarry Smith           y_t = ly[(n10 % (m*n))/m];
107347c6ae99SBarry Smith           /* z_t = z; */
107447c6ae99SBarry Smith           s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
107547c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
107647c6ae99SBarry Smith         }
107747c6ae99SBarry Smith         if (n11 >= 0) { /* right below */
107847c6ae99SBarry Smith           x_t = lx[n11 % m]*dof;
107947c6ae99SBarry Smith           y_t = ly[(n11 % (m*n))/m];
108047c6ae99SBarry Smith           /* z_t = z; */
108147c6ae99SBarry Smith           s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
108247c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
108347c6ae99SBarry Smith         }
108447c6ae99SBarry Smith       }
108547c6ae99SBarry Smith 
108647c6ae99SBarry Smith       for (i=0; i<y; i++) {
108747c6ae99SBarry Smith         if (n12 >= 0) { /* directly left */
108847c6ae99SBarry Smith           x_t = lx[n12 % m]*dof;
108947c6ae99SBarry Smith           y_t = y;
109047c6ae99SBarry Smith           /* z_t = z; */
109147c6ae99SBarry Smith           s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
109247c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
109347c6ae99SBarry Smith         }
109447c6ae99SBarry Smith 
109547c6ae99SBarry Smith         /* Interior */
109647c6ae99SBarry Smith         s_t = bases[rank] + i*x + k*x*y;
109747c6ae99SBarry Smith         for (j=0; j<x; j++) { idx[nn++] = s_t++;}
109847c6ae99SBarry Smith 
109947c6ae99SBarry Smith         if (n14 >= 0) { /* directly right */
110047c6ae99SBarry Smith           x_t = lx[n14 % m]*dof;
110147c6ae99SBarry Smith           y_t = y;
110247c6ae99SBarry Smith           /* z_t = z; */
110347c6ae99SBarry Smith           s_t = bases[n14] + i*x_t + k*x_t*y_t;
110447c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
110547c6ae99SBarry Smith         }
110647c6ae99SBarry Smith       }
110747c6ae99SBarry Smith 
110847c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
110947c6ae99SBarry Smith         if (n15 >= 0) { /* left above */
111047c6ae99SBarry Smith           x_t = lx[n15 % m]*dof;
111147c6ae99SBarry Smith           y_t = ly[(n15 % (m*n))/m];
111247c6ae99SBarry Smith           /* z_t = z; */
111347c6ae99SBarry Smith           s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
111447c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
111547c6ae99SBarry Smith         }
111647c6ae99SBarry Smith         if (n16 >= 0) { /* directly above */
111747c6ae99SBarry Smith           x_t = x;
111847c6ae99SBarry Smith           y_t = ly[(n16 % (m*n))/m];
111947c6ae99SBarry Smith           /* z_t = z; */
112047c6ae99SBarry Smith           s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
112147c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
112247c6ae99SBarry Smith         }
112347c6ae99SBarry Smith         if (n17 >= 0) { /* right above */
112447c6ae99SBarry Smith           x_t = lx[n17 % m]*dof;
112547c6ae99SBarry Smith           y_t = ly[(n17 % (m*n))/m];
112647c6ae99SBarry Smith           /* z_t = z; */
112747c6ae99SBarry Smith           s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
112847c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
112947c6ae99SBarry Smith         }
113047c6ae99SBarry Smith       }
113147c6ae99SBarry Smith     }
113247c6ae99SBarry Smith 
113347c6ae99SBarry Smith     /* Upper Level */
113447c6ae99SBarry Smith     for (k=0; k<s_z; k++) {
113547c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
113647c6ae99SBarry Smith         if (n18 >= 0) { /* left below */
113747c6ae99SBarry Smith           x_t = lx[n18 % m]*dof;
113847c6ae99SBarry Smith           y_t = ly[(n18 % (m*n))/m];
113947c6ae99SBarry Smith           /* z_t = lz[n18 / (m*n)]; */
114047c6ae99SBarry Smith           s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
114147c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
114247c6ae99SBarry Smith         }
114347c6ae99SBarry Smith         if (n19 >= 0) { /* directly below */
114447c6ae99SBarry Smith           x_t = x;
114547c6ae99SBarry Smith           y_t = ly[(n19 % (m*n))/m];
114647c6ae99SBarry Smith           /* z_t = lz[n19 / (m*n)]; */
114747c6ae99SBarry Smith           s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
114847c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
114947c6ae99SBarry Smith         }
115047c6ae99SBarry Smith         if (n20 >= 0) { /* right below */
115147c6ae99SBarry Smith           x_t = lx[n20 % m]*dof;
115247c6ae99SBarry Smith           y_t = ly[(n20 % (m*n))/m];
115347c6ae99SBarry Smith           /* z_t = lz[n20 / (m*n)]; */
115447c6ae99SBarry Smith           s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
115547c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
115647c6ae99SBarry Smith         }
115747c6ae99SBarry Smith       }
115847c6ae99SBarry Smith 
115947c6ae99SBarry Smith       for (i=0; i<y; i++) {
116047c6ae99SBarry Smith         if (n21 >= 0) { /* directly left */
116147c6ae99SBarry Smith           x_t = lx[n21 % m]*dof;
116247c6ae99SBarry Smith           y_t = y;
116347c6ae99SBarry Smith           /* z_t = lz[n21 / (m*n)]; */
116447c6ae99SBarry Smith           s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
116547c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
116647c6ae99SBarry Smith         }
116747c6ae99SBarry Smith 
116847c6ae99SBarry Smith         if (n22 >= 0) { /* middle */
116947c6ae99SBarry Smith           x_t = x;
117047c6ae99SBarry Smith           y_t = y;
117147c6ae99SBarry Smith           /* z_t = lz[n22 / (m*n)]; */
117247c6ae99SBarry Smith           s_t = bases[n22] + i*x_t + k*x_t*y_t;
117347c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
117447c6ae99SBarry Smith         }
117547c6ae99SBarry Smith 
117647c6ae99SBarry Smith         if (n23 >= 0) { /* directly right */
117747c6ae99SBarry Smith           x_t = lx[n23 % m]*dof;
117847c6ae99SBarry Smith           y_t = y;
117947c6ae99SBarry Smith           /* z_t = lz[n23 / (m*n)]; */
118047c6ae99SBarry Smith           s_t = bases[n23] + i*x_t + k*x_t*y_t;
118147c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
118247c6ae99SBarry Smith         }
118347c6ae99SBarry Smith       }
118447c6ae99SBarry Smith 
118547c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
118647c6ae99SBarry Smith         if (n24 >= 0) { /* left above */
118747c6ae99SBarry Smith           x_t = lx[n24 % m]*dof;
118847c6ae99SBarry Smith           y_t = ly[(n24 % (m*n))/m];
118947c6ae99SBarry Smith           /* z_t = lz[n24 / (m*n)]; */
119047c6ae99SBarry Smith           s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
119147c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
119247c6ae99SBarry Smith         }
119347c6ae99SBarry Smith         if (n25 >= 0) { /* directly above */
119447c6ae99SBarry Smith           x_t = x;
119547c6ae99SBarry Smith           y_t = ly[(n25 % (m*n))/m];
119647c6ae99SBarry Smith           /* z_t = lz[n25 / (m*n)]; */
119747c6ae99SBarry Smith           s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
119847c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
119947c6ae99SBarry Smith         }
120047c6ae99SBarry Smith         if (n26 >= 0) { /* right above */
120147c6ae99SBarry Smith           x_t = lx[n26 % m]*dof;
120247c6ae99SBarry Smith           y_t = ly[(n26 % (m*n))/m];
120347c6ae99SBarry Smith           /* z_t = lz[n26 / (m*n)]; */
120447c6ae99SBarry Smith           s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
120547c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
120647c6ae99SBarry Smith         }
120747c6ae99SBarry Smith       }
120847c6ae99SBarry Smith     }
120947c6ae99SBarry Smith   }
121047c6ae99SBarry Smith   /* redo idx to include "missing" ghost points */
121147c6ae99SBarry Smith   /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
121247c6ae99SBarry Smith 
121347c6ae99SBarry Smith   /* Assume Nodes are Internal to the Cube */
121447c6ae99SBarry Smith 
121547c6ae99SBarry Smith   n0  = rank - m*n - m - 1;
121647c6ae99SBarry Smith   n1  = rank - m*n - m;
121747c6ae99SBarry Smith   n2  = rank - m*n - m + 1;
121847c6ae99SBarry Smith   n3  = rank - m*n -1;
121947c6ae99SBarry Smith   n4  = rank - m*n;
122047c6ae99SBarry Smith   n5  = rank - m*n + 1;
122147c6ae99SBarry Smith   n6  = rank - m*n + m - 1;
122247c6ae99SBarry Smith   n7  = rank - m*n + m;
122347c6ae99SBarry Smith   n8  = rank - m*n + m + 1;
122447c6ae99SBarry Smith 
122547c6ae99SBarry Smith   n9  = rank - m - 1;
122647c6ae99SBarry Smith   n10 = rank - m;
122747c6ae99SBarry Smith   n11 = rank - m + 1;
122847c6ae99SBarry Smith   n12 = rank - 1;
122947c6ae99SBarry Smith   n14 = rank + 1;
123047c6ae99SBarry Smith   n15 = rank + m - 1;
123147c6ae99SBarry Smith   n16 = rank + m;
123247c6ae99SBarry Smith   n17 = rank + m + 1;
123347c6ae99SBarry Smith 
123447c6ae99SBarry Smith   n18 = rank + m*n - m - 1;
123547c6ae99SBarry Smith   n19 = rank + m*n - m;
123647c6ae99SBarry Smith   n20 = rank + m*n - m + 1;
123747c6ae99SBarry Smith   n21 = rank + m*n - 1;
123847c6ae99SBarry Smith   n22 = rank + m*n;
123947c6ae99SBarry Smith   n23 = rank + m*n + 1;
124047c6ae99SBarry Smith   n24 = rank + m*n + m - 1;
124147c6ae99SBarry Smith   n25 = rank + m*n + m;
124247c6ae99SBarry Smith   n26 = rank + m*n + m + 1;
124347c6ae99SBarry Smith 
124447c6ae99SBarry Smith   /* Assume Pieces are on Faces of Cube */
124547c6ae99SBarry Smith 
124647c6ae99SBarry Smith   if (xs == 0) { /* First assume not corner or edge */
124747c6ae99SBarry Smith     n0  = rank       -1 - (m*n);
124847c6ae99SBarry Smith     n3  = rank + m   -1 - (m*n);
124947c6ae99SBarry Smith     n6  = rank + 2*m -1 - (m*n);
125047c6ae99SBarry Smith     n9  = rank       -1;
125147c6ae99SBarry Smith     n12 = rank + m   -1;
125247c6ae99SBarry Smith     n15 = rank + 2*m -1;
125347c6ae99SBarry Smith     n18 = rank       -1 + (m*n);
125447c6ae99SBarry Smith     n21 = rank + m   -1 + (m*n);
125547c6ae99SBarry Smith     n24 = rank + 2*m -1 + (m*n);
125647c6ae99SBarry Smith    }
125747c6ae99SBarry Smith 
125847c6ae99SBarry Smith   if (xe == M*dof) { /* First assume not corner or edge */
125947c6ae99SBarry Smith     n2  = rank -2*m +1 - (m*n);
126047c6ae99SBarry Smith     n5  = rank - m  +1 - (m*n);
126147c6ae99SBarry Smith     n8  = rank      +1 - (m*n);
126247c6ae99SBarry Smith     n11 = rank -2*m +1;
126347c6ae99SBarry Smith     n14 = rank - m  +1;
126447c6ae99SBarry Smith     n17 = rank      +1;
126547c6ae99SBarry Smith     n20 = rank -2*m +1 + (m*n);
126647c6ae99SBarry Smith     n23 = rank - m  +1 + (m*n);
126747c6ae99SBarry Smith     n26 = rank      +1 + (m*n);
126847c6ae99SBarry Smith   }
126947c6ae99SBarry Smith 
127047c6ae99SBarry Smith   if (ys==0) { /* First assume not corner or edge */
127147c6ae99SBarry Smith     n0  = rank + m * (n-1) -1 - (m*n);
127247c6ae99SBarry Smith     n1  = rank + m * (n-1)    - (m*n);
127347c6ae99SBarry Smith     n2  = rank + m * (n-1) +1 - (m*n);
127447c6ae99SBarry Smith     n9  = rank + m * (n-1) -1;
127547c6ae99SBarry Smith     n10 = rank + m * (n-1);
127647c6ae99SBarry Smith     n11 = rank + m * (n-1) +1;
127747c6ae99SBarry Smith     n18 = rank + m * (n-1) -1 + (m*n);
127847c6ae99SBarry Smith     n19 = rank + m * (n-1)    + (m*n);
127947c6ae99SBarry Smith     n20 = rank + m * (n-1) +1 + (m*n);
128047c6ae99SBarry Smith   }
128147c6ae99SBarry Smith 
128247c6ae99SBarry Smith   if (ye == N) { /* First assume not corner or edge */
128347c6ae99SBarry Smith     n6  = rank - m * (n-1) -1 - (m*n);
128447c6ae99SBarry Smith     n7  = rank - m * (n-1)    - (m*n);
128547c6ae99SBarry Smith     n8  = rank - m * (n-1) +1 - (m*n);
128647c6ae99SBarry Smith     n15 = rank - m * (n-1) -1;
128747c6ae99SBarry Smith     n16 = rank - m * (n-1);
128847c6ae99SBarry Smith     n17 = rank - m * (n-1) +1;
128947c6ae99SBarry Smith     n24 = rank - m * (n-1) -1 + (m*n);
129047c6ae99SBarry Smith     n25 = rank - m * (n-1)    + (m*n);
129147c6ae99SBarry Smith     n26 = rank - m * (n-1) +1 + (m*n);
129247c6ae99SBarry Smith   }
129347c6ae99SBarry Smith 
129447c6ae99SBarry Smith   if (zs == 0) { /* First assume not corner or edge */
129547c6ae99SBarry Smith     n0 = size - (m*n) + rank - m - 1;
129647c6ae99SBarry Smith     n1 = size - (m*n) + rank - m;
129747c6ae99SBarry Smith     n2 = size - (m*n) + rank - m + 1;
129847c6ae99SBarry Smith     n3 = size - (m*n) + rank - 1;
129947c6ae99SBarry Smith     n4 = size - (m*n) + rank;
130047c6ae99SBarry Smith     n5 = size - (m*n) + rank + 1;
130147c6ae99SBarry Smith     n6 = size - (m*n) + rank + m - 1;
130247c6ae99SBarry Smith     n7 = size - (m*n) + rank + m ;
130347c6ae99SBarry Smith     n8 = size - (m*n) + rank + m + 1;
130447c6ae99SBarry Smith   }
130547c6ae99SBarry Smith 
130647c6ae99SBarry Smith   if (ze == P) { /* First assume not corner or edge */
130747c6ae99SBarry Smith     n18 = (m*n) - (size-rank) - m - 1;
130847c6ae99SBarry Smith     n19 = (m*n) - (size-rank) - m;
130947c6ae99SBarry Smith     n20 = (m*n) - (size-rank) - m + 1;
131047c6ae99SBarry Smith     n21 = (m*n) - (size-rank) - 1;
131147c6ae99SBarry Smith     n22 = (m*n) - (size-rank);
131247c6ae99SBarry Smith     n23 = (m*n) - (size-rank) + 1;
131347c6ae99SBarry Smith     n24 = (m*n) - (size-rank) + m - 1;
131447c6ae99SBarry Smith     n25 = (m*n) - (size-rank) + m;
131547c6ae99SBarry Smith     n26 = (m*n) - (size-rank) + m + 1;
131647c6ae99SBarry Smith   }
131747c6ae99SBarry Smith 
131847c6ae99SBarry Smith   if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */
131947c6ae99SBarry Smith     n0 = size - m*n + rank + m-1 - m;
132047c6ae99SBarry Smith     n3 = size - m*n + rank + m-1;
132147c6ae99SBarry Smith     n6 = size - m*n + rank + m-1 + m;
132247c6ae99SBarry Smith   }
132347c6ae99SBarry Smith 
132447c6ae99SBarry Smith   if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */
132547c6ae99SBarry Smith     n18 = m*n - (size - rank) + m-1 - m;
132647c6ae99SBarry Smith     n21 = m*n - (size - rank) + m-1;
132747c6ae99SBarry Smith     n24 = m*n - (size - rank) + m-1 + m;
132847c6ae99SBarry Smith   }
132947c6ae99SBarry Smith 
133047c6ae99SBarry Smith   if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */
133147c6ae99SBarry Smith     n0  = rank + m*n -1 - m*n;
133247c6ae99SBarry Smith     n9  = rank + m*n -1;
133347c6ae99SBarry Smith     n18 = rank + m*n -1 + m*n;
133447c6ae99SBarry Smith   }
133547c6ae99SBarry Smith 
133647c6ae99SBarry Smith   if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */
133747c6ae99SBarry Smith     n6  = rank - m*(n-1) + m-1 - m*n;
133847c6ae99SBarry Smith     n15 = rank - m*(n-1) + m-1;
133947c6ae99SBarry Smith     n24 = rank - m*(n-1) + m-1 + m*n;
134047c6ae99SBarry Smith   }
134147c6ae99SBarry Smith 
134247c6ae99SBarry Smith   if ((xe==M*dof) && (zs==0)) { /* Assume an edge, not corner */
134347c6ae99SBarry Smith     n2 = size - (m*n-rank) - (m-1) - m;
134447c6ae99SBarry Smith     n5 = size - (m*n-rank) - (m-1);
134547c6ae99SBarry Smith     n8 = size - (m*n-rank) - (m-1) + m;
134647c6ae99SBarry Smith   }
134747c6ae99SBarry Smith 
134847c6ae99SBarry Smith   if ((xe==M*dof) && (ze==P)) { /* Assume an edge, not corner */
134947c6ae99SBarry Smith     n20 = m*n - (size - rank) - (m-1) - m;
135047c6ae99SBarry Smith     n23 = m*n - (size - rank) - (m-1);
135147c6ae99SBarry Smith     n26 = m*n - (size - rank) - (m-1) + m;
135247c6ae99SBarry Smith   }
135347c6ae99SBarry Smith 
135447c6ae99SBarry Smith   if ((xe==M*dof) && (ys==0)) { /* Assume an edge, not corner */
135547c6ae99SBarry Smith     n2  = rank + m*(n-1) - (m-1) - m*n;
135647c6ae99SBarry Smith     n11 = rank + m*(n-1) - (m-1);
135747c6ae99SBarry Smith     n20 = rank + m*(n-1) - (m-1) + m*n;
135847c6ae99SBarry Smith   }
135947c6ae99SBarry Smith 
136047c6ae99SBarry Smith   if ((xe==M*dof) && (ye==N)) { /* Assume an edge, not corner */
136147c6ae99SBarry Smith     n8  = rank - m*n +1 - m*n;
136247c6ae99SBarry Smith     n17 = rank - m*n +1;
136347c6ae99SBarry Smith     n26 = rank - m*n +1 + m*n;
136447c6ae99SBarry Smith   }
136547c6ae99SBarry Smith 
136647c6ae99SBarry Smith   if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */
136747c6ae99SBarry Smith     n0 = size - m + rank -1;
136847c6ae99SBarry Smith     n1 = size - m + rank;
136947c6ae99SBarry Smith     n2 = size - m + rank +1;
137047c6ae99SBarry Smith   }
137147c6ae99SBarry Smith 
137247c6ae99SBarry Smith   if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */
137347c6ae99SBarry Smith     n18 = m*n - (size - rank) + m*(n-1) -1;
137447c6ae99SBarry Smith     n19 = m*n - (size - rank) + m*(n-1);
137547c6ae99SBarry Smith     n20 = m*n - (size - rank) + m*(n-1) +1;
137647c6ae99SBarry Smith   }
137747c6ae99SBarry Smith 
137847c6ae99SBarry Smith   if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */
137947c6ae99SBarry Smith     n6 = size - (m*n-rank) - m * (n-1) -1;
138047c6ae99SBarry Smith     n7 = size - (m*n-rank) - m * (n-1);
138147c6ae99SBarry Smith     n8 = size - (m*n-rank) - m * (n-1) +1;
138247c6ae99SBarry Smith   }
138347c6ae99SBarry Smith 
138447c6ae99SBarry Smith   if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */
138547c6ae99SBarry Smith     n24 = rank - (size-m) -1;
138647c6ae99SBarry Smith     n25 = rank - (size-m);
138747c6ae99SBarry Smith     n26 = rank - (size-m) +1;
138847c6ae99SBarry Smith   }
138947c6ae99SBarry Smith 
139047c6ae99SBarry Smith   /* Check for Corners */
139147c6ae99SBarry Smith   if ((xs==0)   && (ys==0) && (zs==0)) { n0  = size -1;}
139247c6ae99SBarry Smith   if ((xs==0)   && (ys==0) && (ze==P)) { n18 = m*n-1;}
139347c6ae99SBarry Smith   if ((xs==0)   && (ye==N) && (zs==0)) { n6  = (size-1)-m*(n-1);}
139447c6ae99SBarry Smith   if ((xs==0)   && (ye==N) && (ze==P)) { n24 = m-1;}
139547c6ae99SBarry Smith   if ((xe==M*dof) && (ys==0) && (zs==0)) { n2  = size-m;}
139647c6ae99SBarry Smith   if ((xe==M*dof) && (ys==0) && (ze==P)) { n20 = m*n-m;}
139747c6ae99SBarry Smith   if ((xe==M*dof) && (ye==N) && (zs==0)) { n8  = size-m*n;}
139847c6ae99SBarry Smith   if ((xe==M*dof) && (ye==N) && (ze==P)) { n26 = 0;}
139947c6ae99SBarry Smith 
140047c6ae99SBarry Smith   /* Check for when not X,Y, and Z Periodic */
140147c6ae99SBarry Smith 
140247c6ae99SBarry Smith   /* If not X periodic */
140347c6ae99SBarry Smith   if (!DAXPeriodic(wrap)){
140447c6ae99SBarry Smith     if (xs==0)   {n0  = n3  = n6  = n9  = n12 = n15 = n18 = n21 = n24 = -2;}
140547c6ae99SBarry Smith     if (xe==M*dof) {n2  = n5  = n8  = n11 = n14 = n17 = n20 = n23 = n26 = -2;}
140647c6ae99SBarry Smith   }
140747c6ae99SBarry Smith 
140847c6ae99SBarry Smith   /* If not Y periodic */
140947c6ae99SBarry Smith   if (!DAYPeriodic(wrap)){
141047c6ae99SBarry Smith     if (ys==0)   {n0  = n1  = n2  = n9  = n10 = n11 = n18 = n19 = n20 = -2;}
141147c6ae99SBarry Smith     if (ye==N)   {n6  = n7  = n8  = n15 = n16 = n17 = n24 = n25 = n26 = -2;}
141247c6ae99SBarry Smith   }
141347c6ae99SBarry Smith 
141447c6ae99SBarry Smith   /* If not Z periodic */
141547c6ae99SBarry Smith   if (!DAZPeriodic(wrap)){
141647c6ae99SBarry Smith     if (zs==0)   {n0  = n1  = n2  = n3  = n4  = n5  = n6  = n7  = n8  = -2;}
141747c6ae99SBarry Smith     if (ze==P)   {n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;}
141847c6ae99SBarry Smith   }
141947c6ae99SBarry Smith 
142047c6ae99SBarry Smith   nn = 0;
142147c6ae99SBarry Smith 
142247c6ae99SBarry Smith   /* Bottom Level */
142347c6ae99SBarry Smith   for (k=0; k<s_z; k++) {
142447c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
142547c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
142647c6ae99SBarry Smith         x_t = lx[n0 % m]*dof;
142747c6ae99SBarry Smith         y_t = ly[(n0 % (m*n))/m];
142847c6ae99SBarry Smith         z_t = lz[n0 / (m*n)];
142947c6ae99SBarry Smith         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;
143047c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
143147c6ae99SBarry Smith       }
143247c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
143347c6ae99SBarry Smith         x_t = x;
143447c6ae99SBarry Smith         y_t = ly[(n1 % (m*n))/m];
143547c6ae99SBarry Smith         z_t = lz[n1 / (m*n)];
143647c6ae99SBarry Smith         s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
143747c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
143847c6ae99SBarry Smith       }
143947c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
144047c6ae99SBarry Smith         x_t = lx[n2 % m]*dof;
144147c6ae99SBarry Smith         y_t = ly[(n2 % (m*n))/m];
144247c6ae99SBarry Smith         z_t = lz[n2 / (m*n)];
144347c6ae99SBarry Smith         s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
144447c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
144547c6ae99SBarry Smith       }
144647c6ae99SBarry Smith     }
144747c6ae99SBarry Smith 
144847c6ae99SBarry Smith     for (i=0; i<y; i++) {
144947c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
145047c6ae99SBarry Smith         x_t = lx[n3 % m]*dof;
145147c6ae99SBarry Smith         y_t = y;
145247c6ae99SBarry Smith         z_t = lz[n3 / (m*n)];
145347c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
145447c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
145547c6ae99SBarry Smith       }
145647c6ae99SBarry Smith 
145747c6ae99SBarry Smith       if (n4 >= 0) { /* middle */
145847c6ae99SBarry Smith         x_t = x;
145947c6ae99SBarry Smith         y_t = y;
146047c6ae99SBarry Smith         z_t = lz[n4 / (m*n)];
146147c6ae99SBarry Smith         s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
146247c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
146347c6ae99SBarry Smith       }
146447c6ae99SBarry Smith 
146547c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
146647c6ae99SBarry Smith         x_t = lx[n5 % m]*dof;
146747c6ae99SBarry Smith         y_t = y;
146847c6ae99SBarry Smith         z_t = lz[n5 / (m*n)];
146947c6ae99SBarry Smith         s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
147047c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
147147c6ae99SBarry Smith       }
147247c6ae99SBarry Smith     }
147347c6ae99SBarry Smith 
147447c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
147547c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
147647c6ae99SBarry Smith         x_t = lx[n6 % m]*dof;
147747c6ae99SBarry Smith         y_t = ly[(n6 % (m*n))/m];
147847c6ae99SBarry Smith         z_t = lz[n6 / (m*n)];
147947c6ae99SBarry Smith         s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
148047c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
148147c6ae99SBarry Smith       }
148247c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
148347c6ae99SBarry Smith         x_t = x;
148447c6ae99SBarry Smith         y_t = ly[(n7 % (m*n))/m];
148547c6ae99SBarry Smith         z_t = lz[n7 / (m*n)];
148647c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
148747c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
148847c6ae99SBarry Smith       }
148947c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
149047c6ae99SBarry Smith         x_t = lx[n8 % m]*dof;
149147c6ae99SBarry Smith         y_t = ly[(n8 % (m*n))/m];
149247c6ae99SBarry Smith         z_t = lz[n8 / (m*n)];
149347c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
149447c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
149547c6ae99SBarry Smith       }
149647c6ae99SBarry Smith     }
149747c6ae99SBarry Smith   }
149847c6ae99SBarry Smith 
149947c6ae99SBarry Smith   /* Middle Level */
150047c6ae99SBarry Smith   for (k=0; k<z; k++) {
150147c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
150247c6ae99SBarry Smith       if (n9 >= 0) { /* left below */
150347c6ae99SBarry Smith         x_t = lx[n9 % m]*dof;
150447c6ae99SBarry Smith         y_t = ly[(n9 % (m*n))/m];
150547c6ae99SBarry Smith         /* z_t = z; */
150647c6ae99SBarry Smith         s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
150747c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
150847c6ae99SBarry Smith       }
150947c6ae99SBarry Smith       if (n10 >= 0) { /* directly below */
151047c6ae99SBarry Smith         x_t = x;
151147c6ae99SBarry Smith         y_t = ly[(n10 % (m*n))/m];
151247c6ae99SBarry Smith         /* z_t = z; */
151347c6ae99SBarry Smith         s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
151447c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
151547c6ae99SBarry Smith       }
151647c6ae99SBarry Smith       if (n11 >= 0) { /* right below */
151747c6ae99SBarry Smith         x_t = lx[n11 % m]*dof;
151847c6ae99SBarry Smith         y_t = ly[(n11 % (m*n))/m];
151947c6ae99SBarry Smith         /* z_t = z; */
152047c6ae99SBarry Smith         s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
152147c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
152247c6ae99SBarry Smith       }
152347c6ae99SBarry Smith     }
152447c6ae99SBarry Smith 
152547c6ae99SBarry Smith     for (i=0; i<y; i++) {
152647c6ae99SBarry Smith       if (n12 >= 0) { /* directly left */
152747c6ae99SBarry Smith         x_t = lx[n12 % m]*dof;
152847c6ae99SBarry Smith         y_t = y;
152947c6ae99SBarry Smith         /* z_t = z; */
153047c6ae99SBarry Smith         s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
153147c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
153247c6ae99SBarry Smith       }
153347c6ae99SBarry Smith 
153447c6ae99SBarry Smith       /* Interior */
153547c6ae99SBarry Smith       s_t = bases[rank] + i*x + k*x*y;
153647c6ae99SBarry Smith       for (j=0; j<x; j++) { idx[nn++] = s_t++;}
153747c6ae99SBarry Smith 
153847c6ae99SBarry Smith       if (n14 >= 0) { /* directly right */
153947c6ae99SBarry Smith         x_t = lx[n14 % m]*dof;
154047c6ae99SBarry Smith         y_t = y;
154147c6ae99SBarry Smith         /* z_t = z; */
154247c6ae99SBarry Smith         s_t = bases[n14] + i*x_t + k*x_t*y_t;
154347c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
154447c6ae99SBarry Smith       }
154547c6ae99SBarry Smith     }
154647c6ae99SBarry Smith 
154747c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
154847c6ae99SBarry Smith       if (n15 >= 0) { /* left above */
154947c6ae99SBarry Smith         x_t = lx[n15 % m]*dof;
155047c6ae99SBarry Smith         y_t = ly[(n15 % (m*n))/m];
155147c6ae99SBarry Smith         /* z_t = z; */
155247c6ae99SBarry Smith         s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
155347c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
155447c6ae99SBarry Smith       }
155547c6ae99SBarry Smith       if (n16 >= 0) { /* directly above */
155647c6ae99SBarry Smith         x_t = x;
155747c6ae99SBarry Smith         y_t = ly[(n16 % (m*n))/m];
155847c6ae99SBarry Smith         /* z_t = z; */
155947c6ae99SBarry Smith         s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
156047c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
156147c6ae99SBarry Smith       }
156247c6ae99SBarry Smith       if (n17 >= 0) { /* right above */
156347c6ae99SBarry Smith         x_t = lx[n17 % m]*dof;
156447c6ae99SBarry Smith         y_t = ly[(n17 % (m*n))/m];
156547c6ae99SBarry Smith         /* z_t = z; */
156647c6ae99SBarry Smith         s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
156747c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
156847c6ae99SBarry Smith       }
156947c6ae99SBarry Smith     }
157047c6ae99SBarry Smith   }
157147c6ae99SBarry Smith 
157247c6ae99SBarry Smith   /* Upper Level */
157347c6ae99SBarry Smith   for (k=0; k<s_z; k++) {
157447c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
157547c6ae99SBarry Smith       if (n18 >= 0) { /* left below */
157647c6ae99SBarry Smith         x_t = lx[n18 % m]*dof;
157747c6ae99SBarry Smith         y_t = ly[(n18 % (m*n))/m];
157847c6ae99SBarry Smith         /* z_t = lz[n18 / (m*n)]; */
157947c6ae99SBarry Smith         s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
158047c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
158147c6ae99SBarry Smith       }
158247c6ae99SBarry Smith       if (n19 >= 0) { /* directly below */
158347c6ae99SBarry Smith         x_t = x;
158447c6ae99SBarry Smith         y_t = ly[(n19 % (m*n))/m];
158547c6ae99SBarry Smith         /* z_t = lz[n19 / (m*n)]; */
158647c6ae99SBarry Smith         s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
158747c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
158847c6ae99SBarry Smith       }
158947c6ae99SBarry Smith       if (n20 >= 0) { /* right belodof */
159047c6ae99SBarry Smith         x_t = lx[n20 % m]*dof;
159147c6ae99SBarry Smith         y_t = ly[(n20 % (m*n))/m];
159247c6ae99SBarry Smith         /* z_t = lz[n20 / (m*n)]; */
159347c6ae99SBarry Smith         s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
159447c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
159547c6ae99SBarry Smith       }
159647c6ae99SBarry Smith     }
159747c6ae99SBarry Smith 
159847c6ae99SBarry Smith     for (i=0; i<y; i++) {
159947c6ae99SBarry Smith       if (n21 >= 0) { /* directly left */
160047c6ae99SBarry Smith         x_t = lx[n21 % m]*dof;
160147c6ae99SBarry Smith         y_t = y;
160247c6ae99SBarry Smith         /* z_t = lz[n21 / (m*n)]; */
160347c6ae99SBarry Smith         s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
160447c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
160547c6ae99SBarry Smith       }
160647c6ae99SBarry Smith 
160747c6ae99SBarry Smith       if (n22 >= 0) { /* middle */
160847c6ae99SBarry Smith         x_t = x;
160947c6ae99SBarry Smith         y_t = y;
161047c6ae99SBarry Smith         /* z_t = lz[n22 / (m*n)]; */
161147c6ae99SBarry Smith         s_t = bases[n22] + i*x_t + k*x_t*y_t;
161247c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
161347c6ae99SBarry Smith       }
161447c6ae99SBarry Smith 
161547c6ae99SBarry Smith       if (n23 >= 0) { /* directly right */
161647c6ae99SBarry Smith         x_t = lx[n23 % m]*dof;
161747c6ae99SBarry Smith         y_t = y;
161847c6ae99SBarry Smith         /* z_t = lz[n23 / (m*n)]; */
161947c6ae99SBarry Smith         s_t = bases[n23] + i*x_t + k*x_t*y_t;
162047c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
162147c6ae99SBarry Smith       }
162247c6ae99SBarry Smith     }
162347c6ae99SBarry Smith 
162447c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
162547c6ae99SBarry Smith       if (n24 >= 0) { /* left above */
162647c6ae99SBarry Smith         x_t = lx[n24 % m]*dof;
162747c6ae99SBarry Smith         y_t = ly[(n24 % (m*n))/m];
162847c6ae99SBarry Smith         /* z_t = lz[n24 / (m*n)]; */
162947c6ae99SBarry Smith         s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
163047c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
163147c6ae99SBarry Smith       }
163247c6ae99SBarry Smith       if (n25 >= 0) { /* directly above */
163347c6ae99SBarry Smith         x_t = x;
163447c6ae99SBarry Smith         y_t = ly[(n25 % (m*n))/m];
163547c6ae99SBarry Smith         /* z_t = lz[n25 / (m*n)]; */
163647c6ae99SBarry Smith         s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
163747c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
163847c6ae99SBarry Smith       }
163947c6ae99SBarry Smith       if (n26 >= 0) { /* right above */
164047c6ae99SBarry Smith         x_t = lx[n26 % m]*dof;
164147c6ae99SBarry Smith         y_t = ly[(n26 % (m*n))/m];
164247c6ae99SBarry Smith         /* z_t = lz[n26 / (m*n)]; */
164347c6ae99SBarry Smith         s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
164447c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
164547c6ae99SBarry Smith       }
164647c6ae99SBarry Smith     }
164747c6ae99SBarry Smith   }
164847c6ae99SBarry Smith   ierr = PetscFree2(bases,ldims);CHKERRQ(ierr);
164947c6ae99SBarry Smith   dd->gtol      = gtol;
165047c6ae99SBarry Smith   dd->ltog      = ltog;
165147c6ae99SBarry Smith   dd->idx       = idx;
165247c6ae99SBarry Smith   dd->Nl        = nn;
165347c6ae99SBarry Smith   dd->base      = base;
16549a42bb27SBarry Smith   da->ops->view = DMView_DA_3d;
165547c6ae99SBarry Smith 
165647c6ae99SBarry Smith   /*
165747c6ae99SBarry Smith      Set the local to global ordering in the global vector, this allows use
165847c6ae99SBarry Smith      of VecSetValuesLocal().
165947c6ae99SBarry Smith   */
166047c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,nn,idx,PETSC_OWN_POINTER,&dd->ltogmap);CHKERRQ(ierr);
166147c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingBlock(dd->ltogmap,dd->w,&dd->ltogmapb);CHKERRQ(ierr);
166247c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,dd->ltogmap);CHKERRQ(ierr);
166347c6ae99SBarry Smith 
166447c6ae99SBarry Smith   dd->ltol = PETSC_NULL;
166547c6ae99SBarry Smith   dd->ao   = PETSC_NULL;
166647c6ae99SBarry Smith   PetscFunctionReturn(0);
166747c6ae99SBarry Smith }
1668*564755cdSBarry Smith 
166947c6ae99SBarry Smith 
167047c6ae99SBarry Smith #undef __FUNCT__
167147c6ae99SBarry Smith #define __FUNCT__ "DACreate3d"
167247c6ae99SBarry Smith /*@C
167347c6ae99SBarry Smith    DACreate3d - Creates an object that will manage the communication of three-dimensional
167447c6ae99SBarry Smith    regular array data that is distributed across some processors.
167547c6ae99SBarry Smith 
167647c6ae99SBarry Smith    Collective on MPI_Comm
167747c6ae99SBarry Smith 
167847c6ae99SBarry Smith    Input Parameters:
167947c6ae99SBarry Smith +  comm - MPI communicator
168047c6ae99SBarry Smith .  wrap - type of periodicity the array should have, if any.  Use one
168147c6ae99SBarry Smith           of DA_NONPERIODIC, DA_XPERIODIC, DA_YPERIODIC, DA_ZPERIODIC, DA_XYPERIODIC, DA_XZPERIODIC, DA_YZPERIODIC, DA_XYZPERIODIC, or DA_XYZGHOSTED.
168247c6ae99SBarry Smith .  stencil_type - Type of stencil (DA_STENCIL_STAR or DA_STENCIL_BOX)
168347c6ae99SBarry Smith .  M,N,P - global dimension in each direction of the array (use -M, -N, and or -P to indicate that it may be set to a different value
168447c6ae99SBarry Smith             from the command line with -da_grid_x <M> -da_grid_y <N> -da_grid_z <P>)
168547c6ae99SBarry Smith .  m,n,p - corresponding number of processors in each dimension
168647c6ae99SBarry Smith            (or PETSC_DECIDE to have calculated)
168747c6ae99SBarry Smith .  dof - number of degrees of freedom per node
168847c6ae99SBarry Smith .  lx, ly, lz - arrays containing the number of nodes in each cell along
168947c6ae99SBarry Smith           the x, y, and z coordinates, or PETSC_NULL. If non-null, these
169047c6ae99SBarry Smith           must be of length as m,n,p and the corresponding
169147c6ae99SBarry Smith           m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of
169247c6ae99SBarry Smith           the ly[] must N, sum of the lz[] must be P
169347c6ae99SBarry Smith -  s - stencil width
169447c6ae99SBarry Smith 
169547c6ae99SBarry Smith    Output Parameter:
169647c6ae99SBarry Smith .  da - the resulting distributed array object
169747c6ae99SBarry Smith 
169847c6ae99SBarry Smith    Options Database Key:
16999a42bb27SBarry Smith +  -da_view - Calls DMView() at the conclusion of DACreate3d()
170047c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
170147c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
170247c6ae99SBarry Smith .  -da_grid_z <nz> - number of grid points in z direction, if P < 0
170347c6ae99SBarry Smith .  -da_processors_x <MX> - number of processors in x direction
170447c6ae99SBarry Smith .  -da_processors_y <MY> - number of processors in y direction
170547c6ae99SBarry Smith .  -da_processors_z <MZ> - number of processors in z direction
170647c6ae99SBarry Smith .  -da_refine_x - refinement ratio in x direction
170747c6ae99SBarry Smith .  -da_refine_y - refinement ratio in y direction
170847c6ae99SBarry Smith -  -da_refine_y - refinement ratio in z direction
170947c6ae99SBarry Smith 
171047c6ae99SBarry Smith    Level: beginner
171147c6ae99SBarry Smith 
171247c6ae99SBarry Smith    Notes:
171347c6ae99SBarry Smith    The stencil type DA_STENCIL_STAR with width 1 corresponds to the
171447c6ae99SBarry Smith    standard 7-pt stencil, while DA_STENCIL_BOX with width 1 denotes
171547c6ae99SBarry Smith    the standard 27-pt stencil.
171647c6ae99SBarry Smith 
171747c6ae99SBarry Smith    The array data itself is NOT stored in the DA, it is stored in Vec objects;
1718*564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
1719*564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
172047c6ae99SBarry Smith 
172147c6ae99SBarry Smith .keywords: distributed array, create, three-dimensional
172247c6ae99SBarry Smith 
17239a42bb27SBarry Smith .seealso: DMDestroy(), DMView(), DACreate1d(), DACreate2d(), DMGlobalToLocalBegin(), DAGetRefinementFactor(),
17249a42bb27SBarry Smith           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DALocalToLocalBegin(), DALocalToLocalEnd(), DASetRefinementFactor(),
1725*564755cdSBarry Smith           DAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DACreateNaturalVector(), DALoad(), DAGetOwnershipRanges()
172647c6ae99SBarry Smith 
172747c6ae99SBarry Smith @*/
172847c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DACreate3d(MPI_Comm comm,DAPeriodicType wrap,DAStencilType stencil_type,PetscInt M,
17299a42bb27SBarry Smith                PetscInt N,PetscInt P,PetscInt m,PetscInt n,PetscInt p,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],const PetscInt lz[],DM *da)
173047c6ae99SBarry Smith {
173147c6ae99SBarry Smith   PetscErrorCode ierr;
173247c6ae99SBarry Smith 
173347c6ae99SBarry Smith   PetscFunctionBegin;
173447c6ae99SBarry Smith   ierr = DACreate(comm, da);CHKERRQ(ierr);
173547c6ae99SBarry Smith   ierr = DASetDim(*da, 3);CHKERRQ(ierr);
173647c6ae99SBarry Smith   ierr = DASetSizes(*da, M, N, P);CHKERRQ(ierr);
173747c6ae99SBarry Smith   ierr = DASetNumProcs(*da, m, n, p);CHKERRQ(ierr);
173847c6ae99SBarry Smith   ierr = DASetPeriodicity(*da, wrap);CHKERRQ(ierr);
173947c6ae99SBarry Smith   ierr = DASetDof(*da, dof);CHKERRQ(ierr);
174047c6ae99SBarry Smith   ierr = DASetStencilType(*da, stencil_type);CHKERRQ(ierr);
174147c6ae99SBarry Smith   ierr = DASetStencilWidth(*da, s);CHKERRQ(ierr);
174247c6ae99SBarry Smith   ierr = DASetOwnershipRanges(*da, lx, ly, lz);CHKERRQ(ierr);
174347c6ae99SBarry Smith   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
17449a42bb27SBarry Smith   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
17459a42bb27SBarry Smith   ierr = DMSetUp(*da);CHKERRQ(ierr);
174647c6ae99SBarry Smith   PetscFunctionReturn(0);
174747c6ae99SBarry Smith }
1748