xref: /petsc/src/dm/impls/da/da3.c (revision 8135c375e96a7bf70d67c02f02cd826154817056)
17d0a6c19SBarry Smith 
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 
7af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h>     /*I   "petscdmda.h"    I*/
847c6ae99SBarry Smith 
99804daf3SBarry Smith #include <petscdraw.h>
10e0877f53SBarry Smith static PetscErrorCode DMView_DA_3d(DM da,PetscViewer viewer)
1147c6ae99SBarry Smith {
1247c6ae99SBarry Smith   PetscErrorCode ierr;
1347c6ae99SBarry Smith   PetscMPIInt    rank;
149a42bb27SBarry Smith   PetscBool      iascii,isdraw,isbinary;
1547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
169a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
179a42bb27SBarry Smith   PetscBool ismatlab;
189a42bb27SBarry Smith #endif
1947c6ae99SBarry Smith 
2047c6ae99SBarry Smith   PetscFunctionBegin;
21ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr);
2247c6ae99SBarry Smith 
23251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
24251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
25251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
269a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
27251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr);
289a42bb27SBarry Smith #endif
2947c6ae99SBarry Smith   if (iascii) {
3047c6ae99SBarry Smith     PetscViewerFormat format;
3147c6ae99SBarry Smith 
321575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
3347c6ae99SBarry Smith     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
34*8135c375SStefano Zampini     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL && format != PETSC_VIEWER_ASCII_GLVIS) {
35aa219208SBarry Smith       DMDALocalInfo info;
36aa219208SBarry Smith       ierr = DMDAGetLocalInfo(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)
416636e97aSMatthew G Knepley       if (da->coordinates) {
4247c6ae99SBarry Smith         PetscInt        last;
4347c6ae99SBarry Smith         const PetscReal *coors;
446636e97aSMatthew G Knepley         ierr = VecGetArrayRead(da->coordinates,&coors);CHKERRQ(ierr);
456636e97aSMatthew G Knepley         ierr = VecGetLocalSize(da->coordinates,&last);CHKERRQ(ierr);
4647c6ae99SBarry Smith         last = last - 3;
4757622a8eSBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Lower left corner %g %g %g : Upper right %g %g %g\n",(double)coors[0],(double)coors[1],(double)coors[2],(double)coors[last],(double)coors[last+1],(double)coors[last+2]);CHKERRQ(ierr);
486636e97aSMatthew G Knepley         ierr = VecRestoreArrayRead(da->coordinates,&coors);CHKERRQ(ierr);
4947c6ae99SBarry Smith       }
5047c6ae99SBarry Smith #endif
5147c6ae99SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
521575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
53*8135c375SStefano Zampini     } else if (format == PETSC_VIEWER_ASCII_GLVIS) {
54*8135c375SStefano Zampini       ierr = DMView_DA_GLVis(da,viewer);CHKERRQ(ierr);
55616157d6SJed Brown     } else {
56616157d6SJed Brown       ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr);
5747c6ae99SBarry Smith     }
5847c6ae99SBarry Smith   } else if (isdraw) {
5947c6ae99SBarry Smith     PetscDraw      draw;
6047c6ae99SBarry Smith     PetscReal      ymin = -1.0,ymax = (PetscReal)dd->N;
6147c6ae99SBarry Smith     PetscReal      xmin = -1.0,xmax = (PetscReal)((dd->M+2)*dd->P),x,y,ycoord,xcoord;
628ea3bf28SBarry Smith     PetscInt       k,plane,base;
638ea3bf28SBarry Smith     const PetscInt *idx;
6447c6ae99SBarry Smith     char           node[10];
6547c6ae99SBarry Smith     PetscBool      isnull;
6647c6ae99SBarry Smith 
6747c6ae99SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
685b399a63SLisandro Dalcin     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
695b399a63SLisandro Dalcin     if (isnull) PetscFunctionReturn(0);
7047c6ae99SBarry Smith 
715b399a63SLisandro Dalcin     ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
725b399a63SLisandro Dalcin     ierr = PetscDrawClear(draw);CHKERRQ(ierr);
735b399a63SLisandro Dalcin     ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
745b399a63SLisandro Dalcin 
755b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
7647c6ae99SBarry Smith     /* first processor draw all node lines */
7747c6ae99SBarry Smith     if (!rank) {
7847c6ae99SBarry Smith       for (k=0; k<dd->P; k++) {
7947c6ae99SBarry Smith         ymin = 0.0; ymax = (PetscReal)(dd->N - 1);
8047c6ae99SBarry Smith         for (xmin=(PetscReal)(k*(dd->M+1)); xmin<(PetscReal)(dd->M+(k*(dd->M+1))); xmin++) {
8147c6ae99SBarry Smith           ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
8247c6ae99SBarry Smith         }
8347c6ae99SBarry Smith         xmin = (PetscReal)(k*(dd->M+1)); xmax = xmin + (PetscReal)(dd->M - 1);
8447c6ae99SBarry Smith         for (ymin=0; ymin<(PetscReal)dd->N; ymin++) {
8547c6ae99SBarry Smith           ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
8647c6ae99SBarry Smith         }
8747c6ae99SBarry Smith       }
8847c6ae99SBarry Smith     }
895b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
905b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
9147c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
9247c6ae99SBarry Smith 
935b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
945b399a63SLisandro Dalcin     /*Go through and draw for each plane*/
955b399a63SLisandro Dalcin     for (k=0; k<dd->P; k++) {
9647c6ae99SBarry Smith       if ((k >= dd->zs) && (k < dd->ze)) {
9747c6ae99SBarry Smith         /* draw my box */
9847c6ae99SBarry Smith         ymin = dd->ys;
9947c6ae99SBarry Smith         ymax = dd->ye - 1;
10047c6ae99SBarry Smith         xmin = dd->xs/dd->w    + (dd->M+1)*k;
10147c6ae99SBarry Smith         xmax =(dd->xe-1)/dd->w + (dd->M+1)*k;
10247c6ae99SBarry Smith 
10347c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
10447c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
10547c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
10647c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
10747c6ae99SBarry Smith 
10847c6ae99SBarry Smith         xmin = dd->xs/dd->w;
10947c6ae99SBarry Smith         xmax =(dd->xe-1)/dd->w;
11047c6ae99SBarry Smith 
111832b7cebSLisandro Dalcin         /* identify which processor owns the box */
112832b7cebSLisandro Dalcin         ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)rank);CHKERRQ(ierr);
113832b7cebSLisandro Dalcin         ierr = PetscDrawString(draw,xmin+(dd->M+1)*k+.2,ymin+.3,PETSC_DRAW_RED,node);CHKERRQ(ierr);
11447c6ae99SBarry Smith         /* put in numbers*/
11547c6ae99SBarry Smith         base = (dd->base+(dd->xe-dd->xs)*(dd->ye-dd->ys)*(k-dd->zs))/dd->w;
11647c6ae99SBarry Smith         for (y=ymin; y<=ymax; y++) {
11747c6ae99SBarry Smith           for (x=xmin+(dd->M+1)*k; x<=xmax+(dd->M+1)*k; x++) {
118832b7cebSLisandro Dalcin             ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)base++);CHKERRQ(ierr);
11947c6ae99SBarry Smith             ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
12047c6ae99SBarry Smith           }
12147c6ae99SBarry Smith         }
12247c6ae99SBarry Smith 
12347c6ae99SBarry Smith       }
12447c6ae99SBarry Smith     }
1255b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
1265b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
12747c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
12847c6ae99SBarry Smith 
1295b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
13047c6ae99SBarry Smith     for (k=0-dd->s; k<dd->P+dd->s; k++) {
13147c6ae99SBarry Smith       /* Go through and draw for each plane */
13247c6ae99SBarry Smith       if ((k >= dd->Zs) && (k < dd->Ze)) {
13347c6ae99SBarry Smith         /* overlay ghost numbers, useful for error checking */
134565245c5SBarry Smith         base = (dd->Xe-dd->Xs)*(dd->Ye-dd->Ys)*(k-dd->Zs)/dd->w;
13545b6f7e9SBarry Smith         ierr = ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr);
13647c6ae99SBarry Smith         plane=k;
137565245c5SBarry Smith         /* Keep z wrap around points on the drawing */
1388865f1eaSKarl Rupp         if (k<0) plane=dd->P+k;
1398865f1eaSKarl Rupp         if (k>=dd->P) plane=k-dd->P;
14047c6ae99SBarry Smith         ymin = dd->Ys; ymax = dd->Ye;
14147c6ae99SBarry Smith         xmin = (dd->M+1)*plane*dd->w;
14247c6ae99SBarry Smith         xmax = (dd->M+1)*plane*dd->w+dd->M*dd->w;
14347c6ae99SBarry Smith         for (y=ymin; y<ymax; y++) {
14447c6ae99SBarry Smith           for (x=xmin+dd->Xs; x<xmin+dd->Xe; x+=dd->w) {
1458ea3bf28SBarry Smith             sprintf(node,"%d",(int)(idx[base]));
14647c6ae99SBarry Smith             ycoord = y;
14747c6ae99SBarry Smith             /*Keep y wrap around points on drawing */
1488865f1eaSKarl Rupp             if (y<0) ycoord = dd->N+y;
1498865f1eaSKarl Rupp             if (y>=dd->N) ycoord = y-dd->N;
15047c6ae99SBarry Smith             xcoord = x;   /* Keep x wrap points on drawing */
1518865f1eaSKarl Rupp             if (x<xmin) xcoord = xmax - (xmin-x);
1528865f1eaSKarl Rupp             if (x>=xmax) xcoord = xmin + (x-xmax);
15347c6ae99SBarry Smith             ierr = PetscDrawString(draw,xcoord/dd->w,ycoord,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
154565245c5SBarry Smith             base++;
15547c6ae99SBarry Smith           }
15647c6ae99SBarry Smith         }
15745b6f7e9SBarry Smith         ierr = ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr);
15847c6ae99SBarry Smith       }
15947c6ae99SBarry Smith     }
1605b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
1615b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
16247c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
163832b7cebSLisandro Dalcin     ierr = PetscDrawSave(draw);CHKERRQ(ierr);
1649a42bb27SBarry Smith   } else if (isbinary) {
1659a42bb27SBarry Smith     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
1669a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1679a42bb27SBarry Smith   } else if (ismatlab) {
1689a42bb27SBarry Smith     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
1699a42bb27SBarry Smith #endif
17011aeaf0aSBarry Smith   }
17147c6ae99SBarry Smith   PetscFunctionReturn(0);
17247c6ae99SBarry Smith }
17347c6ae99SBarry Smith 
1747087cfbeSBarry Smith PetscErrorCode  DMSetUp_DA_3D(DM da)
17547c6ae99SBarry Smith {
17647c6ae99SBarry Smith   DM_DA            *dd          = (DM_DA*)da->data;
17747c6ae99SBarry Smith   const PetscInt   M            = dd->M;
17847c6ae99SBarry Smith   const PetscInt   N            = dd->N;
17947c6ae99SBarry Smith   const PetscInt   P            = dd->P;
18047c6ae99SBarry Smith   PetscInt         m            = dd->m;
18147c6ae99SBarry Smith   PetscInt         n            = dd->n;
18247c6ae99SBarry Smith   PetscInt         p            = dd->p;
18347c6ae99SBarry Smith   const PetscInt   dof          = dd->w;
18447c6ae99SBarry Smith   const PetscInt   s            = dd->s;
185bff4a2f0SMatthew G. Knepley   DMBoundaryType   bx           = dd->bx;
186bff4a2f0SMatthew G. Knepley   DMBoundaryType   by           = dd->by;
187bff4a2f0SMatthew G. Knepley   DMBoundaryType   bz           = dd->bz;
18819fd82e9SBarry Smith   DMDAStencilType  stencil_type = dd->stencil_type;
18947c6ae99SBarry Smith   PetscInt         *lx          = dd->lx;
19047c6ae99SBarry Smith   PetscInt         *ly          = dd->ly;
19147c6ae99SBarry Smith   PetscInt         *lz          = dd->lz;
19247c6ae99SBarry Smith   MPI_Comm         comm;
19347c6ae99SBarry Smith   PetscMPIInt      rank,size;
194ce00eea3SSatish Balay   PetscInt         xs = 0,xe,ys = 0,ye,zs = 0,ze,x = 0,y = 0,z = 0;
195bd1fc5aeSBarry Smith   PetscInt         Xs,Xe,Ys,Ye,Zs,Ze,IXs,IXe,IYs,IYe,IZs,IZe,pm;
1968ea3bf28SBarry Smith   PetscInt         left,right,up,down,bottom,top,i,j,k,*idx,nn;
19747c6ae99SBarry Smith   PetscInt         n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14;
19847c6ae99SBarry Smith   PetscInt         n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26;
199db87c5ecSEthan Coon   PetscInt         *bases,*ldims,base,x_t,y_t,z_t,s_t,count,s_x,s_y,s_z;
20047c6ae99SBarry Smith   PetscInt         sn0  = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0;
20147c6ae99SBarry Smith   PetscInt         sn8  = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0;
20247c6ae99SBarry Smith   PetscInt         sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0;
20347c6ae99SBarry Smith   Vec              local,global;
204bd1fc5aeSBarry Smith   VecScatter       gtol;
20545b6f7e9SBarry Smith   IS               to,from;
2066f951b95Secoon   PetscBool        twod;
20747c6ae99SBarry Smith   PetscErrorCode   ierr;
20847c6ae99SBarry Smith 
2096f951b95Secoon 
21047c6ae99SBarry Smith   PetscFunctionBegin;
211bff4a2f0SMatthew G. Knepley   if (stencil_type == DMDA_STENCIL_BOX && (bx == DM_BOUNDARY_MIRROR || by == DM_BOUNDARY_MIRROR || bz == DM_BOUNDARY_MIRROR)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Mirror boundary and box stencil");
21240088605SBarry Smith   if ((bx == DM_BOUNDARY_MIRROR) || (by == DM_BOUNDARY_MIRROR) || (bz == DM_BOUNDARY_MIRROR)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Mirror boundary not supported yet in 3d");
21347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr);
2143855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES)
215bafee8b4SSatish Balay   if (((PetscInt64) M)*((PetscInt64) N)*((PetscInt64) P)*((PetscInt64) dof) > (PetscInt64) PETSC_MPI_INT_MAX) SETERRQ3(comm,PETSC_ERR_INT_OVERFLOW,"Mesh of %D by %D by %D (dof) is too large for 32 bit indices",M,N,dof);
2163855c12bSBarry Smith #endif
2173855c12bSBarry Smith 
21847c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
21947c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
22047c6ae99SBarry Smith 
22147c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
22247c6ae99SBarry Smith     if (m < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
22347c6ae99SBarry Smith     else if (m > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
22447c6ae99SBarry Smith   }
22547c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
22647c6ae99SBarry Smith     if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
22747c6ae99SBarry Smith     else if (n > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
22847c6ae99SBarry Smith   }
22947c6ae99SBarry Smith   if (p != PETSC_DECIDE) {
23047c6ae99SBarry Smith     if (p < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %D",p);
23147c6ae99SBarry Smith     else if (p > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %D %d",p,size);
23247c6ae99SBarry Smith   }
2330716a85fSBarry Smith   if ((m > 0) && (n > 0) && (p > 0) && (m*n*p != size)) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m %D * n %D * p %D != size %d",m,n,p,size);
23447c6ae99SBarry Smith 
23547c6ae99SBarry Smith   /* Partition the array among the processors */
23647c6ae99SBarry Smith   if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) {
23747c6ae99SBarry Smith     m = size/(n*p);
23847c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
23947c6ae99SBarry Smith     n = size/(m*p);
24047c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
24147c6ae99SBarry Smith     p = size/(m*n);
24247c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
24347c6ae99SBarry Smith     /* try for squarish distribution */
244369cc0aeSBarry Smith     m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N*p)));
24547c6ae99SBarry Smith     if (!m) m = 1;
24647c6ae99SBarry Smith     while (m > 0) {
24747c6ae99SBarry Smith       n = size/(m*p);
24847c6ae99SBarry Smith       if (m*n*p == size) break;
24947c6ae99SBarry Smith       m--;
25047c6ae99SBarry Smith     }
25147c6ae99SBarry Smith     if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %D",p);
25247c6ae99SBarry Smith     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
25347c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
25447c6ae99SBarry Smith     /* try for squarish distribution */
255369cc0aeSBarry Smith     m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n)));
25647c6ae99SBarry Smith     if (!m) m = 1;
25747c6ae99SBarry Smith     while (m > 0) {
25847c6ae99SBarry Smith       p = size/(m*n);
25947c6ae99SBarry Smith       if (m*n*p == size) break;
26047c6ae99SBarry Smith       m--;
26147c6ae99SBarry Smith     }
26247c6ae99SBarry Smith     if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %D",n);
26347c6ae99SBarry Smith     if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;}
26447c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
26547c6ae99SBarry Smith     /* try for squarish distribution */
266369cc0aeSBarry Smith     n = (int)(0.5 + PetscSqrtReal(((PetscReal)N)*((PetscReal)size)/((PetscReal)P*m)));
26747c6ae99SBarry Smith     if (!n) n = 1;
26847c6ae99SBarry Smith     while (n > 0) {
26947c6ae99SBarry Smith       p = size/(m*n);
27047c6ae99SBarry Smith       if (m*n*p == size) break;
27147c6ae99SBarry Smith       n--;
27247c6ae99SBarry Smith     }
27347c6ae99SBarry Smith     if (!n) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %D",n);
27447c6ae99SBarry Smith     if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;}
27547c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
27647c6ae99SBarry Smith     /* try for squarish distribution */
277369cc0aeSBarry Smith     n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),(PetscReal)(1./3.)));
27847c6ae99SBarry Smith     if (!n) n = 1;
27947c6ae99SBarry Smith     while (n > 0) {
28047c6ae99SBarry Smith       pm = size/n;
28147c6ae99SBarry Smith       if (n*pm == size) break;
28247c6ae99SBarry Smith       n--;
28347c6ae99SBarry Smith     }
28447c6ae99SBarry Smith     if (!n) n = 1;
285369cc0aeSBarry Smith     m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n)));
28647c6ae99SBarry Smith     if (!m) m = 1;
28747c6ae99SBarry Smith     while (m > 0) {
28847c6ae99SBarry Smith       p = size/(m*n);
28947c6ae99SBarry Smith       if (m*n*p == size) break;
29047c6ae99SBarry Smith       m--;
29147c6ae99SBarry Smith     }
29247c6ae99SBarry Smith     if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;}
293ce94432eSBarry Smith   } else if (m*n*p != size) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
29447c6ae99SBarry Smith 
295ce94432eSBarry Smith   if (m*n*p != size) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_PLIB,"Could not find good partition");
296ce94432eSBarry Smith   if (M < m) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
297ce94432eSBarry Smith   if (N < n) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
298ce94432eSBarry Smith   if (P < p) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %D %D",P,p);
29947c6ae99SBarry Smith 
30047c6ae99SBarry Smith   /*
30147c6ae99SBarry Smith      Determine locally owned region
30247c6ae99SBarry Smith      [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes
30347c6ae99SBarry Smith   */
30447c6ae99SBarry Smith 
30547c6ae99SBarry Smith   if (!lx) {
306785e854fSJed Brown     ierr = PetscMalloc1(m, &dd->lx);CHKERRQ(ierr);
30747c6ae99SBarry Smith     lx   = dd->lx;
3088865f1eaSKarl Rupp     for (i=0; i<m; i++) lx[i] = M/m + ((M % m) > (i % m));
30947c6ae99SBarry Smith   }
31047c6ae99SBarry Smith   x  = lx[rank % m];
31147c6ae99SBarry Smith   xs = 0;
3128865f1eaSKarl Rupp   for (i=0; i<(rank%m); i++) xs += lx[i];
313bff4a2f0SMatthew G. Knepley   if ((x < s) && ((m > 1) || (bx == DM_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s);
31447c6ae99SBarry Smith 
31547c6ae99SBarry Smith   if (!ly) {
316785e854fSJed Brown     ierr = PetscMalloc1(n, &dd->ly);CHKERRQ(ierr);
31747c6ae99SBarry Smith     ly   = dd->ly;
3188865f1eaSKarl Rupp     for (i=0; i<n; i++) ly[i] = N/n + ((N % n) > (i % n));
31947c6ae99SBarry Smith   }
32047c6ae99SBarry Smith   y = ly[(rank % (m*n))/m];
321bff4a2f0SMatthew G. Knepley   if ((y < s) && ((n > 1) || (by == DM_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s);
32230729d88SBarry Smith 
32347c6ae99SBarry Smith   ys = 0;
3248865f1eaSKarl Rupp   for (i=0; i<(rank % (m*n))/m; i++) ys += ly[i];
32547c6ae99SBarry Smith 
32647c6ae99SBarry Smith   if (!lz) {
327785e854fSJed Brown     ierr = PetscMalloc1(p, &dd->lz);CHKERRQ(ierr);
32847c6ae99SBarry Smith     lz = dd->lz;
3298865f1eaSKarl Rupp     for (i=0; i<p; i++) lz[i] = P/p + ((P % p) > (i % p));
33047c6ae99SBarry Smith   }
33147c6ae99SBarry Smith   z = lz[rank/(m*n)];
332bcea557cSEthan Coon 
333fdc81ce8SEthan Coon   /* note this is different than x- and y-, as we will handle as an important special
334bff4a2f0SMatthew G. Knepley    case when p=P=1 and DM_BOUNDARY_PERIODIC and s > z.  This is to deal with 2D problems
335fdc81ce8SEthan Coon    in a 3D code.  Additional code for this case is noted with "2d case" comments */
3366f951b95Secoon   twod = PETSC_FALSE;
3378865f1eaSKarl Rupp   if (P == 1) twod = PETSC_TRUE;
338bff4a2f0SMatthew G. Knepley   else if ((z < s) && ((p > 1) || (bz == DM_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local z-width of domain z %D is smaller than stencil width s %D",z,s);
33947c6ae99SBarry Smith   zs = 0;
3408865f1eaSKarl Rupp   for (i=0; i<(rank/(m*n)); i++) zs += lz[i];
34147c6ae99SBarry Smith   ye = ys + y;
34247c6ae99SBarry Smith   xe = xs + x;
34347c6ae99SBarry Smith   ze = zs + z;
34447c6ae99SBarry Smith 
34588661749SPeter Brune   /* determine ghost region (Xs) and region scattered into (IXs)  */
346d9c9ebe5SPeter Brune   if (xs-s > 0) {
347d9c9ebe5SPeter Brune     Xs = xs - s; IXs = xs - s;
34888661749SPeter Brune   } else {
3498865f1eaSKarl Rupp     if (bx) Xs = xs - s;
3508865f1eaSKarl Rupp     else Xs = 0;
35188661749SPeter Brune     IXs = 0;
35288661749SPeter Brune   }
353d9c9ebe5SPeter Brune   if (xe+s <= M) {
354d9c9ebe5SPeter Brune     Xe = xe + s; IXe = xe + s;
35588661749SPeter Brune   } else {
35688661749SPeter Brune     if (bx) {
357d9c9ebe5SPeter Brune       Xs = xs - s; Xe = xe + s;
3588865f1eaSKarl Rupp     } else Xe = M;
35988661749SPeter Brune     IXe = M;
36088661749SPeter Brune   }
36147c6ae99SBarry Smith 
362bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) {
363d9c9ebe5SPeter Brune     IXs = xs - s;
364d9c9ebe5SPeter Brune     IXe = xe + s;
365d9c9ebe5SPeter Brune     Xs  = xs - s;
366d9c9ebe5SPeter Brune     Xe  = xe + s;
36788661749SPeter Brune   }
36888661749SPeter Brune 
369d9c9ebe5SPeter Brune   if (ys-s > 0) {
370d9c9ebe5SPeter Brune     Ys = ys - s; IYs = ys - s;
37188661749SPeter Brune   } else {
3728865f1eaSKarl Rupp     if (by) Ys = ys - s;
3738865f1eaSKarl Rupp     else Ys = 0;
37488661749SPeter Brune     IYs = 0;
37588661749SPeter Brune   }
376d9c9ebe5SPeter Brune   if (ye+s <= N) {
377d9c9ebe5SPeter Brune     Ye = ye + s; IYe = ye + s;
37888661749SPeter Brune   } else {
3798865f1eaSKarl Rupp     if (by) Ye = ye + s;
3808865f1eaSKarl Rupp     else Ye = N;
38188661749SPeter Brune     IYe = N;
38288661749SPeter Brune   }
38388661749SPeter Brune 
384bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) {
385d9c9ebe5SPeter Brune     IYs = ys - s;
386d9c9ebe5SPeter Brune     IYe = ye + s;
387d9c9ebe5SPeter Brune     Ys  = ys - s;
388d9c9ebe5SPeter Brune     Ye  = ye + s;
38988661749SPeter Brune   }
39088661749SPeter Brune 
391d9c9ebe5SPeter Brune   if (zs-s > 0) {
392d9c9ebe5SPeter Brune     Zs = zs - s; IZs = zs - s;
39388661749SPeter Brune   } else {
3948865f1eaSKarl Rupp     if (bz) Zs = zs - s;
3958865f1eaSKarl Rupp     else Zs = 0;
39688661749SPeter Brune     IZs = 0;
39788661749SPeter Brune   }
398d9c9ebe5SPeter Brune   if (ze+s <= P) {
399d9c9ebe5SPeter Brune     Ze = ze + s; IZe = ze + s;
40088661749SPeter Brune   } else {
4018865f1eaSKarl Rupp     if (bz) Ze = ze + s;
4028865f1eaSKarl Rupp     else Ze = P;
40388661749SPeter Brune     IZe = P;
40488661749SPeter Brune   }
40588661749SPeter Brune 
406bff4a2f0SMatthew G. Knepley   if (bz == DM_BOUNDARY_PERIODIC || bz == DM_BOUNDARY_MIRROR) {
407d9c9ebe5SPeter Brune     IZs = zs - s;
408d9c9ebe5SPeter Brune     IZe = ze + s;
409d9c9ebe5SPeter Brune     Zs  = zs - s;
410d9c9ebe5SPeter Brune     Ze  = ze + s;
41188661749SPeter Brune   }
41247c6ae99SBarry Smith 
41347c6ae99SBarry Smith   /* Resize all X parameters to reflect w */
414d9c9ebe5SPeter Brune   s_x = s;
415d9c9ebe5SPeter Brune   s_y = s;
416d9c9ebe5SPeter Brune   s_z = s;
41747c6ae99SBarry Smith 
41847c6ae99SBarry Smith   /* determine starting point of each processor */
41947c6ae99SBarry Smith   nn       = x*y*z;
420dcca6d9dSJed Brown   ierr     = PetscMalloc2(size+1,&bases,size,&ldims);CHKERRQ(ierr);
42147c6ae99SBarry Smith   ierr     = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
42247c6ae99SBarry Smith   bases[0] = 0;
4238865f1eaSKarl Rupp   for (i=1; i<=size; i++) bases[i] = ldims[i-1];
4248865f1eaSKarl Rupp   for (i=1; i<=size; i++) bases[i] += bases[i-1];
425ce00eea3SSatish Balay   base = bases[rank]*dof;
42647c6ae99SBarry Smith 
42747c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
428ce00eea3SSatish Balay   dd->Nlocal = x*y*z*dof;
429b1fb7eb7SBarry Smith   ierr       = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global);CHKERRQ(ierr);
430ce00eea3SSatish Balay   dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs)*dof;
431b1fb7eb7SBarry Smith   ierr       = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local);CHKERRQ(ierr);
43247c6ae99SBarry Smith 
43347c6ae99SBarry Smith   /* generate appropriate vector scatters */
43447c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
43502fe608eSBarry Smith   ierr   = PetscMalloc1((IXe-IXs)*(IYe-IYs)*(IZe-IZs),&idx);CHKERRQ(ierr);
436ce00eea3SSatish Balay   left   = xs - Xs; right = left + x;
43747c6ae99SBarry Smith   bottom = ys - Ys; top = bottom + y;
43847c6ae99SBarry Smith   down   = zs - Zs; up  = down + z;
43947c6ae99SBarry Smith   count  = 0;
44047c6ae99SBarry Smith   for (i=down; i<up; i++) {
44147c6ae99SBarry Smith     for (j=bottom; j<top; j++) {
442ce00eea3SSatish Balay       for (k=left; k<right; k++) {
443ce00eea3SSatish Balay         idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
44447c6ae99SBarry Smith       }
44547c6ae99SBarry Smith     }
44647c6ae99SBarry Smith   }
44747c6ae99SBarry Smith 
448ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
449ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
450d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_BOX) {
451ce00eea3SSatish Balay     left   = IXs - Xs; right = left + (IXe-IXs);
452ce00eea3SSatish Balay     bottom = IYs - Ys; top = bottom + (IYe-IYs);
453ce00eea3SSatish Balay     down   = IZs - Zs; up  = down + (IZe-IZs);
454ce00eea3SSatish Balay     count  = 0;
455ce00eea3SSatish Balay     for (i=down; i<up; i++) {
456ce00eea3SSatish Balay       for (j=bottom; j<top; j++) {
457ce00eea3SSatish Balay         for (k=left; k<right; k++) {
458ce00eea3SSatish Balay           idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
459ce00eea3SSatish Balay         }
460ce00eea3SSatish Balay       }
461ce00eea3SSatish Balay     }
462ce00eea3SSatish Balay     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
46347c6ae99SBarry Smith   } else {
46447c6ae99SBarry Smith     /* This is way ugly! We need to list the funny cross type region */
465ce00eea3SSatish Balay     left   = xs - Xs; right = left + x;
46647c6ae99SBarry Smith     bottom = ys - Ys; top = bottom + y;
46747c6ae99SBarry Smith     down   = zs - Zs;   up  = down + z;
46847c6ae99SBarry Smith     count  = 0;
469ce00eea3SSatish Balay     /* the bottom chunck */
470ce00eea3SSatish Balay     for (i=(IZs-Zs); i<down; i++) {
47147c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
472ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
47347c6ae99SBarry Smith       }
47447c6ae99SBarry Smith     }
47547c6ae99SBarry Smith     /* the middle piece */
47647c6ae99SBarry Smith     for (i=down; i<up; i++) {
47747c6ae99SBarry Smith       /* front */
478ce00eea3SSatish Balay       for (j=(IYs-Ys); j<bottom; j++) {
479ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
48047c6ae99SBarry Smith       }
48147c6ae99SBarry Smith       /* middle */
48247c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
483ce00eea3SSatish Balay         for (k=IXs-Xs; k<IXe-Xs; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
48447c6ae99SBarry Smith       }
48547c6ae99SBarry Smith       /* back */
486ce00eea3SSatish Balay       for (j=top; j<top+IYe-ye; j++) {
487ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
48847c6ae99SBarry Smith       }
48947c6ae99SBarry Smith     }
49047c6ae99SBarry Smith     /* the top piece */
491ce00eea3SSatish Balay     for (i=up; i<up+IZe-ze; i++) {
49247c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
493ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
49447c6ae99SBarry Smith       }
49547c6ae99SBarry Smith     }
49647c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
49747c6ae99SBarry Smith   }
49847c6ae99SBarry Smith 
49947c6ae99SBarry Smith   /* determine who lies on each side of use stored in    n24 n25 n26
50047c6ae99SBarry Smith                                                          n21 n22 n23
50147c6ae99SBarry Smith                                                          n18 n19 n20
50247c6ae99SBarry Smith 
50347c6ae99SBarry Smith                                                          n15 n16 n17
50447c6ae99SBarry Smith                                                          n12     n14
50547c6ae99SBarry Smith                                                          n9  n10 n11
50647c6ae99SBarry Smith 
50747c6ae99SBarry Smith                                                          n6  n7  n8
50847c6ae99SBarry Smith                                                          n3  n4  n5
50947c6ae99SBarry Smith                                                          n0  n1  n2
51047c6ae99SBarry Smith   */
51147c6ae99SBarry Smith 
51247c6ae99SBarry Smith   /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
51347c6ae99SBarry Smith   /* Assume Nodes are Internal to the Cube */
51447c6ae99SBarry Smith   n0 = rank - m*n - m - 1;
51547c6ae99SBarry Smith   n1 = rank - m*n - m;
51647c6ae99SBarry Smith   n2 = rank - m*n - m + 1;
51747c6ae99SBarry Smith   n3 = rank - m*n -1;
51847c6ae99SBarry Smith   n4 = rank - m*n;
51947c6ae99SBarry Smith   n5 = rank - m*n + 1;
52047c6ae99SBarry Smith   n6 = rank - m*n + m - 1;
52147c6ae99SBarry Smith   n7 = rank - m*n + m;
52247c6ae99SBarry Smith   n8 = rank - m*n + m + 1;
52347c6ae99SBarry Smith 
52447c6ae99SBarry Smith   n9  = rank - m - 1;
52547c6ae99SBarry Smith   n10 = rank - m;
52647c6ae99SBarry Smith   n11 = rank - m + 1;
52747c6ae99SBarry Smith   n12 = rank - 1;
52847c6ae99SBarry Smith   n14 = rank + 1;
52947c6ae99SBarry Smith   n15 = rank + m - 1;
53047c6ae99SBarry Smith   n16 = rank + m;
53147c6ae99SBarry Smith   n17 = rank + m + 1;
53247c6ae99SBarry Smith 
53347c6ae99SBarry Smith   n18 = rank + m*n - m - 1;
53447c6ae99SBarry Smith   n19 = rank + m*n - m;
53547c6ae99SBarry Smith   n20 = rank + m*n - m + 1;
53647c6ae99SBarry Smith   n21 = rank + m*n - 1;
53747c6ae99SBarry Smith   n22 = rank + m*n;
53847c6ae99SBarry Smith   n23 = rank + m*n + 1;
53947c6ae99SBarry Smith   n24 = rank + m*n + m - 1;
54047c6ae99SBarry Smith   n25 = rank + m*n + m;
54147c6ae99SBarry Smith   n26 = rank + m*n + m + 1;
54247c6ae99SBarry Smith 
54347c6ae99SBarry Smith   /* Assume Pieces are on Faces of Cube */
54447c6ae99SBarry Smith 
54547c6ae99SBarry Smith   if (xs == 0) { /* First assume not corner or edge */
54647c6ae99SBarry Smith     n0  = rank       -1 - (m*n);
54747c6ae99SBarry Smith     n3  = rank + m   -1 - (m*n);
54847c6ae99SBarry Smith     n6  = rank + 2*m -1 - (m*n);
54947c6ae99SBarry Smith     n9  = rank       -1;
55047c6ae99SBarry Smith     n12 = rank + m   -1;
55147c6ae99SBarry Smith     n15 = rank + 2*m -1;
55247c6ae99SBarry Smith     n18 = rank       -1 + (m*n);
55347c6ae99SBarry Smith     n21 = rank + m   -1 + (m*n);
55447c6ae99SBarry Smith     n24 = rank + 2*m -1 + (m*n);
55547c6ae99SBarry Smith   }
55647c6ae99SBarry Smith 
557ce00eea3SSatish Balay   if (xe == M) { /* First assume not corner or edge */
55847c6ae99SBarry Smith     n2  = rank -2*m +1 - (m*n);
55947c6ae99SBarry Smith     n5  = rank - m  +1 - (m*n);
56047c6ae99SBarry Smith     n8  = rank      +1 - (m*n);
56147c6ae99SBarry Smith     n11 = rank -2*m +1;
56247c6ae99SBarry Smith     n14 = rank - m  +1;
56347c6ae99SBarry Smith     n17 = rank      +1;
56447c6ae99SBarry Smith     n20 = rank -2*m +1 + (m*n);
56547c6ae99SBarry Smith     n23 = rank - m  +1 + (m*n);
56647c6ae99SBarry Smith     n26 = rank      +1 + (m*n);
56747c6ae99SBarry Smith   }
56847c6ae99SBarry Smith 
56947c6ae99SBarry Smith   if (ys==0) { /* First assume not corner or edge */
57047c6ae99SBarry Smith     n0  = rank + m * (n-1) -1 - (m*n);
57147c6ae99SBarry Smith     n1  = rank + m * (n-1)    - (m*n);
57247c6ae99SBarry Smith     n2  = rank + m * (n-1) +1 - (m*n);
57347c6ae99SBarry Smith     n9  = rank + m * (n-1) -1;
57447c6ae99SBarry Smith     n10 = rank + m * (n-1);
57547c6ae99SBarry Smith     n11 = rank + m * (n-1) +1;
57647c6ae99SBarry Smith     n18 = rank + m * (n-1) -1 + (m*n);
57747c6ae99SBarry Smith     n19 = rank + m * (n-1)    + (m*n);
57847c6ae99SBarry Smith     n20 = rank + m * (n-1) +1 + (m*n);
57947c6ae99SBarry Smith   }
58047c6ae99SBarry Smith 
58147c6ae99SBarry Smith   if (ye == N) { /* First assume not corner or edge */
58247c6ae99SBarry Smith     n6  = rank - m * (n-1) -1 - (m*n);
58347c6ae99SBarry Smith     n7  = rank - m * (n-1)    - (m*n);
58447c6ae99SBarry Smith     n8  = rank - m * (n-1) +1 - (m*n);
58547c6ae99SBarry Smith     n15 = rank - m * (n-1) -1;
58647c6ae99SBarry Smith     n16 = rank - m * (n-1);
58747c6ae99SBarry Smith     n17 = rank - m * (n-1) +1;
58847c6ae99SBarry Smith     n24 = rank - m * (n-1) -1 + (m*n);
58947c6ae99SBarry Smith     n25 = rank - m * (n-1)    + (m*n);
59047c6ae99SBarry Smith     n26 = rank - m * (n-1) +1 + (m*n);
59147c6ae99SBarry Smith   }
59247c6ae99SBarry Smith 
59347c6ae99SBarry Smith   if (zs == 0) { /* First assume not corner or edge */
59447c6ae99SBarry Smith     n0 = size - (m*n) + rank - m - 1;
59547c6ae99SBarry Smith     n1 = size - (m*n) + rank - m;
59647c6ae99SBarry Smith     n2 = size - (m*n) + rank - m + 1;
59747c6ae99SBarry Smith     n3 = size - (m*n) + rank - 1;
59847c6ae99SBarry Smith     n4 = size - (m*n) + rank;
59947c6ae99SBarry Smith     n5 = size - (m*n) + rank + 1;
60047c6ae99SBarry Smith     n6 = size - (m*n) + rank + m - 1;
60147c6ae99SBarry Smith     n7 = size - (m*n) + rank + m;
60247c6ae99SBarry Smith     n8 = size - (m*n) + rank + m + 1;
60347c6ae99SBarry Smith   }
60447c6ae99SBarry Smith 
60547c6ae99SBarry Smith   if (ze == P) { /* First assume not corner or edge */
60647c6ae99SBarry Smith     n18 = (m*n) - (size-rank) - m - 1;
60747c6ae99SBarry Smith     n19 = (m*n) - (size-rank) - m;
60847c6ae99SBarry Smith     n20 = (m*n) - (size-rank) - m + 1;
60947c6ae99SBarry Smith     n21 = (m*n) - (size-rank) - 1;
61047c6ae99SBarry Smith     n22 = (m*n) - (size-rank);
61147c6ae99SBarry Smith     n23 = (m*n) - (size-rank) + 1;
61247c6ae99SBarry Smith     n24 = (m*n) - (size-rank) + m - 1;
61347c6ae99SBarry Smith     n25 = (m*n) - (size-rank) + m;
61447c6ae99SBarry Smith     n26 = (m*n) - (size-rank) + m + 1;
61547c6ae99SBarry Smith   }
61647c6ae99SBarry Smith 
61747c6ae99SBarry Smith   if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */
61847c6ae99SBarry Smith     n0 = size - m*n + rank + m-1 - m;
61947c6ae99SBarry Smith     n3 = size - m*n + rank + m-1;
62047c6ae99SBarry Smith     n6 = size - m*n + rank + m-1 + m;
62147c6ae99SBarry Smith   }
62247c6ae99SBarry Smith 
62347c6ae99SBarry Smith   if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */
62447c6ae99SBarry Smith     n18 = m*n - (size - rank) + m-1 - m;
62547c6ae99SBarry Smith     n21 = m*n - (size - rank) + m-1;
62647c6ae99SBarry Smith     n24 = m*n - (size - rank) + m-1 + m;
62747c6ae99SBarry Smith   }
62847c6ae99SBarry Smith 
62947c6ae99SBarry Smith   if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */
63047c6ae99SBarry Smith     n0  = rank + m*n -1 - m*n;
63147c6ae99SBarry Smith     n9  = rank + m*n -1;
63247c6ae99SBarry Smith     n18 = rank + m*n -1 + m*n;
63347c6ae99SBarry Smith   }
63447c6ae99SBarry Smith 
63547c6ae99SBarry Smith   if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */
63647c6ae99SBarry Smith     n6  = rank - m*(n-1) + m-1 - m*n;
63747c6ae99SBarry Smith     n15 = rank - m*(n-1) + m-1;
63847c6ae99SBarry Smith     n24 = rank - m*(n-1) + m-1 + m*n;
63947c6ae99SBarry Smith   }
64047c6ae99SBarry Smith 
641ce00eea3SSatish Balay   if ((xe==M) && (zs==0)) { /* Assume an edge, not corner */
64247c6ae99SBarry Smith     n2 = size - (m*n-rank) - (m-1) - m;
64347c6ae99SBarry Smith     n5 = size - (m*n-rank) - (m-1);
64447c6ae99SBarry Smith     n8 = size - (m*n-rank) - (m-1) + m;
64547c6ae99SBarry Smith   }
64647c6ae99SBarry Smith 
647ce00eea3SSatish Balay   if ((xe==M) && (ze==P)) { /* Assume an edge, not corner */
64847c6ae99SBarry Smith     n20 = m*n - (size - rank) - (m-1) - m;
64947c6ae99SBarry Smith     n23 = m*n - (size - rank) - (m-1);
65047c6ae99SBarry Smith     n26 = m*n - (size - rank) - (m-1) + m;
65147c6ae99SBarry Smith   }
65247c6ae99SBarry Smith 
653ce00eea3SSatish Balay   if ((xe==M) && (ys==0)) { /* Assume an edge, not corner */
65447c6ae99SBarry Smith     n2  = rank + m*(n-1) - (m-1) - m*n;
65547c6ae99SBarry Smith     n11 = rank + m*(n-1) - (m-1);
65647c6ae99SBarry Smith     n20 = rank + m*(n-1) - (m-1) + m*n;
65747c6ae99SBarry Smith   }
65847c6ae99SBarry Smith 
659ce00eea3SSatish Balay   if ((xe==M) && (ye==N)) { /* Assume an edge, not corner */
66047c6ae99SBarry Smith     n8  = rank - m*n +1 - m*n;
66147c6ae99SBarry Smith     n17 = rank - m*n +1;
66247c6ae99SBarry Smith     n26 = rank - m*n +1 + m*n;
66347c6ae99SBarry Smith   }
66447c6ae99SBarry Smith 
66547c6ae99SBarry Smith   if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */
66647c6ae99SBarry Smith     n0 = size - m + rank -1;
66747c6ae99SBarry Smith     n1 = size - m + rank;
66847c6ae99SBarry Smith     n2 = size - m + rank +1;
66947c6ae99SBarry Smith   }
67047c6ae99SBarry Smith 
67147c6ae99SBarry Smith   if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */
67247c6ae99SBarry Smith     n18 = m*n - (size - rank) + m*(n-1) -1;
67347c6ae99SBarry Smith     n19 = m*n - (size - rank) + m*(n-1);
67447c6ae99SBarry Smith     n20 = m*n - (size - rank) + m*(n-1) +1;
67547c6ae99SBarry Smith   }
67647c6ae99SBarry Smith 
67747c6ae99SBarry Smith   if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */
67847c6ae99SBarry Smith     n6 = size - (m*n-rank) - m * (n-1) -1;
67947c6ae99SBarry Smith     n7 = size - (m*n-rank) - m * (n-1);
68047c6ae99SBarry Smith     n8 = size - (m*n-rank) - m * (n-1) +1;
68147c6ae99SBarry Smith   }
68247c6ae99SBarry Smith 
68347c6ae99SBarry Smith   if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */
68447c6ae99SBarry Smith     n24 = rank - (size-m) -1;
68547c6ae99SBarry Smith     n25 = rank - (size-m);
68647c6ae99SBarry Smith     n26 = rank - (size-m) +1;
68747c6ae99SBarry Smith   }
68847c6ae99SBarry Smith 
68947c6ae99SBarry Smith   /* Check for Corners */
6908865f1eaSKarl Rupp   if ((xs==0) && (ys==0) && (zs==0)) n0  = size -1;
6918865f1eaSKarl Rupp   if ((xs==0) && (ys==0) && (ze==P)) n18 = m*n-1;
6928865f1eaSKarl Rupp   if ((xs==0) && (ye==N) && (zs==0)) n6  = (size-1)-m*(n-1);
6938865f1eaSKarl Rupp   if ((xs==0) && (ye==N) && (ze==P)) n24 = m-1;
6948865f1eaSKarl Rupp   if ((xe==M) && (ys==0) && (zs==0)) n2  = size-m;
6958865f1eaSKarl Rupp   if ((xe==M) && (ys==0) && (ze==P)) n20 = m*n-m;
6968865f1eaSKarl Rupp   if ((xe==M) && (ye==N) && (zs==0)) n8  = size-m*n;
6978865f1eaSKarl Rupp   if ((xe==M) && (ye==N) && (ze==P)) n26 = 0;
69847c6ae99SBarry Smith 
69947c6ae99SBarry Smith   /* Check for when not X,Y, and Z Periodic */
70047c6ae99SBarry Smith 
70147c6ae99SBarry Smith   /* If not X periodic */
702bff4a2f0SMatthew G. Knepley   if (bx != DM_BOUNDARY_PERIODIC) {
7038865f1eaSKarl Rupp     if (xs==0) n0 = n3 = n6 = n9  = n12 = n15 = n18 = n21 = n24 = -2;
7048865f1eaSKarl Rupp     if (xe==M) n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2;
70547c6ae99SBarry Smith   }
70647c6ae99SBarry Smith 
70747c6ae99SBarry Smith   /* If not Y periodic */
708bff4a2f0SMatthew G. Knepley   if (by != DM_BOUNDARY_PERIODIC) {
7098865f1eaSKarl Rupp     if (ys==0) n0 = n1 = n2 = n9  = n10 = n11 = n18 = n19 = n20 = -2;
7108865f1eaSKarl Rupp     if (ye==N) n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2;
71147c6ae99SBarry Smith   }
71247c6ae99SBarry Smith 
71347c6ae99SBarry Smith   /* If not Z periodic */
714bff4a2f0SMatthew G. Knepley   if (bz != DM_BOUNDARY_PERIODIC) {
7158865f1eaSKarl Rupp     if (zs==0) n0  = n1  = n2  = n3  = n4  = n5  = n6  = n7  = n8  = -2;
7168865f1eaSKarl Rupp     if (ze==P) n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;
71747c6ae99SBarry Smith   }
71847c6ae99SBarry Smith 
719785e854fSJed Brown   ierr = PetscMalloc1(27,&dd->neighbors);CHKERRQ(ierr);
7208865f1eaSKarl Rupp 
72147c6ae99SBarry Smith   dd->neighbors[0]  = n0;
72247c6ae99SBarry Smith   dd->neighbors[1]  = n1;
72347c6ae99SBarry Smith   dd->neighbors[2]  = n2;
72447c6ae99SBarry Smith   dd->neighbors[3]  = n3;
72547c6ae99SBarry Smith   dd->neighbors[4]  = n4;
72647c6ae99SBarry Smith   dd->neighbors[5]  = n5;
72747c6ae99SBarry Smith   dd->neighbors[6]  = n6;
72847c6ae99SBarry Smith   dd->neighbors[7]  = n7;
72947c6ae99SBarry Smith   dd->neighbors[8]  = n8;
73047c6ae99SBarry Smith   dd->neighbors[9]  = n9;
73147c6ae99SBarry Smith   dd->neighbors[10] = n10;
73247c6ae99SBarry Smith   dd->neighbors[11] = n11;
73347c6ae99SBarry Smith   dd->neighbors[12] = n12;
73447c6ae99SBarry Smith   dd->neighbors[13] = rank;
73547c6ae99SBarry Smith   dd->neighbors[14] = n14;
73647c6ae99SBarry Smith   dd->neighbors[15] = n15;
73747c6ae99SBarry Smith   dd->neighbors[16] = n16;
73847c6ae99SBarry Smith   dd->neighbors[17] = n17;
73947c6ae99SBarry Smith   dd->neighbors[18] = n18;
74047c6ae99SBarry Smith   dd->neighbors[19] = n19;
74147c6ae99SBarry Smith   dd->neighbors[20] = n20;
74247c6ae99SBarry Smith   dd->neighbors[21] = n21;
74347c6ae99SBarry Smith   dd->neighbors[22] = n22;
74447c6ae99SBarry Smith   dd->neighbors[23] = n23;
74547c6ae99SBarry Smith   dd->neighbors[24] = n24;
74647c6ae99SBarry Smith   dd->neighbors[25] = n25;
74747c6ae99SBarry Smith   dd->neighbors[26] = n26;
74847c6ae99SBarry Smith 
74947c6ae99SBarry Smith   /* If star stencil then delete the corner neighbors */
750d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
75147c6ae99SBarry Smith     /* save information about corner neighbors */
75247c6ae99SBarry Smith     sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7;
75347c6ae99SBarry Smith     sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18;
75447c6ae99SBarry Smith     sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25;
75547c6ae99SBarry Smith     sn26 = n26;
7568865f1eaSKarl Rupp     n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 = n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1;
75747c6ae99SBarry Smith   }
75847c6ae99SBarry Smith 
759785e854fSJed Brown   ierr = PetscMalloc1((Xe-Xs)*(Ye-Ys)*(Ze-Zs),&idx);CHKERRQ(ierr);
76047c6ae99SBarry Smith 
76147c6ae99SBarry Smith   nn = 0;
76247c6ae99SBarry Smith   /* Bottom Level */
76347c6ae99SBarry Smith   for (k=0; k<s_z; k++) {
76447c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
76547c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
766ce00eea3SSatish Balay         x_t = lx[n0 % m];
76747c6ae99SBarry Smith         y_t = ly[(n0 % (m*n))/m];
76847c6ae99SBarry Smith         z_t = lz[n0 / (m*n)];
76947c6ae99SBarry 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;
7708865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n0] + x_t*y_t*z_t - (s_y-i)*x_t - s_x; /* 2D case */
7718865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
77247c6ae99SBarry Smith       }
77347c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
77447c6ae99SBarry Smith         x_t = x;
77547c6ae99SBarry Smith         y_t = ly[(n1 % (m*n))/m];
77647c6ae99SBarry Smith         z_t = lz[n1 / (m*n)];
77747c6ae99SBarry 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;
7788865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */
7798865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
78047c6ae99SBarry Smith       }
78147c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
782ce00eea3SSatish Balay         x_t = lx[n2 % m];
78347c6ae99SBarry Smith         y_t = ly[(n2 % (m*n))/m];
78447c6ae99SBarry Smith         z_t = lz[n2 / (m*n)];
78547c6ae99SBarry 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;
7868865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */
7878865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
78847c6ae99SBarry Smith       }
78947c6ae99SBarry Smith     }
79047c6ae99SBarry Smith 
79147c6ae99SBarry Smith     for (i=0; i<y; i++) {
79247c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
793ce00eea3SSatish Balay         x_t = lx[n3 % m];
79447c6ae99SBarry Smith         y_t = y;
79547c6ae99SBarry Smith         z_t = lz[n3 / (m*n)];
79647c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
7978865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - x_t*y_t; /* 2D case */
7988865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
79947c6ae99SBarry Smith       }
80047c6ae99SBarry Smith 
80147c6ae99SBarry Smith       if (n4 >= 0) { /* middle */
80247c6ae99SBarry Smith         x_t = x;
80347c6ae99SBarry Smith         y_t = y;
80447c6ae99SBarry Smith         z_t = lz[n4 / (m*n)];
80547c6ae99SBarry Smith         s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8068865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n4] + i*x_t + x_t*y_t*z_t - x_t*y_t; /* 2D case */
8078865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
808bff4a2f0SMatthew G. Knepley       } else if (bz == DM_BOUNDARY_MIRROR) {
8098865f1eaSKarl Rupp         for (j=0; j<x; j++) idx[nn++] = 0;
81047c6ae99SBarry Smith       }
81147c6ae99SBarry Smith 
81247c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
813ce00eea3SSatish Balay         x_t = lx[n5 % m];
81447c6ae99SBarry Smith         y_t = y;
81547c6ae99SBarry Smith         z_t = lz[n5 / (m*n)];
81647c6ae99SBarry Smith         s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8178865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n5] + i*x_t + x_t*y_t*z_t - x_t*y_t; /* 2D case */
8188865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
81947c6ae99SBarry Smith       }
82047c6ae99SBarry Smith     }
82147c6ae99SBarry Smith 
82247c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
82347c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
824ce00eea3SSatish Balay         x_t = lx[n6 % m];
82547c6ae99SBarry Smith         y_t = ly[(n6 % (m*n))/m];
82647c6ae99SBarry Smith         z_t = lz[n6 / (m*n)];
82747c6ae99SBarry Smith         s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8288865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - x_t*y_t; /* 2D case */
8298865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
83047c6ae99SBarry Smith       }
83147c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
83247c6ae99SBarry Smith         x_t = x;
83347c6ae99SBarry Smith         y_t = ly[(n7 % (m*n))/m];
83447c6ae99SBarry Smith         z_t = lz[n7 / (m*n)];
83547c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8368865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - x_t*y_t; /* 2D case */
8378865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
83847c6ae99SBarry Smith       }
83947c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
840ce00eea3SSatish Balay         x_t = lx[n8 % m];
84147c6ae99SBarry Smith         y_t = ly[(n8 % (m*n))/m];
84247c6ae99SBarry Smith         z_t = lz[n8 / (m*n)];
84347c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8448865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - x_t*y_t; /* 2D case */
8458865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
84647c6ae99SBarry Smith       }
84747c6ae99SBarry Smith     }
84847c6ae99SBarry Smith   }
84947c6ae99SBarry Smith 
85047c6ae99SBarry Smith   /* Middle Level */
85147c6ae99SBarry Smith   for (k=0; k<z; k++) {
85247c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
85347c6ae99SBarry Smith       if (n9 >= 0) { /* left below */
854ce00eea3SSatish Balay         x_t = lx[n9 % m];
85547c6ae99SBarry Smith         y_t = ly[(n9 % (m*n))/m];
85647c6ae99SBarry Smith         /* z_t = z; */
85747c6ae99SBarry Smith         s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
8588865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
85947c6ae99SBarry Smith       }
86047c6ae99SBarry Smith       if (n10 >= 0) { /* directly below */
86147c6ae99SBarry Smith         x_t = x;
86247c6ae99SBarry Smith         y_t = ly[(n10 % (m*n))/m];
86347c6ae99SBarry Smith         /* z_t = z; */
86447c6ae99SBarry Smith         s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
8658865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
866bff4a2f0SMatthew G. Knepley       }  else if (by == DM_BOUNDARY_MIRROR) {
8678865f1eaSKarl Rupp         for (j=0; j<x; j++) idx[nn++] = 0;
86847c6ae99SBarry Smith       }
86947c6ae99SBarry Smith       if (n11 >= 0) { /* right below */
870ce00eea3SSatish Balay         x_t = lx[n11 % m];
87147c6ae99SBarry Smith         y_t = ly[(n11 % (m*n))/m];
87247c6ae99SBarry Smith         /* z_t = z; */
87347c6ae99SBarry Smith         s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
8748865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
87547c6ae99SBarry Smith       }
87647c6ae99SBarry Smith     }
87747c6ae99SBarry Smith 
87847c6ae99SBarry Smith     for (i=0; i<y; i++) {
87947c6ae99SBarry Smith       if (n12 >= 0) { /* directly left */
880ce00eea3SSatish Balay         x_t = lx[n12 % m];
88147c6ae99SBarry Smith         y_t = y;
88247c6ae99SBarry Smith         /* z_t = z; */
88347c6ae99SBarry Smith         s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
8848865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
885bff4a2f0SMatthew G. Knepley       }  else if (bx == DM_BOUNDARY_MIRROR) {
8868865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = 0;
88747c6ae99SBarry Smith       }
88847c6ae99SBarry Smith 
88947c6ae99SBarry Smith       /* Interior */
89047c6ae99SBarry Smith       s_t = bases[rank] + i*x + k*x*y;
8918865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = s_t++;
89247c6ae99SBarry Smith 
89347c6ae99SBarry Smith       if (n14 >= 0) { /* directly right */
894ce00eea3SSatish Balay         x_t = lx[n14 % m];
89547c6ae99SBarry Smith         y_t = y;
89647c6ae99SBarry Smith         /* z_t = z; */
89747c6ae99SBarry Smith         s_t = bases[n14] + i*x_t + k*x_t*y_t;
8988865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
899bff4a2f0SMatthew G. Knepley       } else if (bx == DM_BOUNDARY_MIRROR) {
9008865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = 0;
90147c6ae99SBarry Smith       }
90247c6ae99SBarry Smith     }
90347c6ae99SBarry Smith 
90447c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
90547c6ae99SBarry Smith       if (n15 >= 0) { /* left above */
906ce00eea3SSatish Balay         x_t = lx[n15 % m];
90747c6ae99SBarry Smith         y_t = ly[(n15 % (m*n))/m];
90847c6ae99SBarry Smith         /* z_t = z; */
90947c6ae99SBarry Smith         s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
9108865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
91147c6ae99SBarry Smith       }
91247c6ae99SBarry Smith       if (n16 >= 0) { /* directly above */
91347c6ae99SBarry Smith         x_t = x;
91447c6ae99SBarry Smith         y_t = ly[(n16 % (m*n))/m];
91547c6ae99SBarry Smith         /* z_t = z; */
91647c6ae99SBarry Smith         s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
9178865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
918bff4a2f0SMatthew G. Knepley       } else if (by == DM_BOUNDARY_MIRROR) {
9198865f1eaSKarl Rupp         for (j=0; j<x; j++) idx[nn++] = 0;
92047c6ae99SBarry Smith       }
92147c6ae99SBarry Smith       if (n17 >= 0) { /* right above */
922ce00eea3SSatish Balay         x_t = lx[n17 % m];
92347c6ae99SBarry Smith         y_t = ly[(n17 % (m*n))/m];
92447c6ae99SBarry Smith         /* z_t = z; */
92547c6ae99SBarry Smith         s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
9268865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
92747c6ae99SBarry Smith       }
92847c6ae99SBarry Smith     }
92947c6ae99SBarry Smith   }
93047c6ae99SBarry Smith 
93147c6ae99SBarry Smith   /* Upper Level */
93247c6ae99SBarry Smith   for (k=0; k<s_z; k++) {
93347c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
93447c6ae99SBarry Smith       if (n18 >= 0) { /* left below */
935ce00eea3SSatish Balay         x_t = lx[n18 % m];
93647c6ae99SBarry Smith         y_t = ly[(n18 % (m*n))/m];
93747c6ae99SBarry Smith         /* z_t = lz[n18 / (m*n)]; */
93847c6ae99SBarry Smith         s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
9398865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n18] - (s_y-i)*x_t -s_x + x_t*y_t; /* 2d case */
9408865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
94147c6ae99SBarry Smith       }
94247c6ae99SBarry Smith       if (n19 >= 0) { /* directly below */
94347c6ae99SBarry Smith         x_t = x;
94447c6ae99SBarry Smith         y_t = ly[(n19 % (m*n))/m];
94547c6ae99SBarry Smith         /* z_t = lz[n19 / (m*n)]; */
94647c6ae99SBarry Smith         s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
9478865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n19] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */
9488865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
94947c6ae99SBarry Smith       }
95047c6ae99SBarry Smith       if (n20 >= 0) { /* right below */
951ce00eea3SSatish Balay         x_t = lx[n20 % m];
95247c6ae99SBarry Smith         y_t = ly[(n20 % (m*n))/m];
95347c6ae99SBarry Smith         /* z_t = lz[n20 / (m*n)]; */
95447c6ae99SBarry Smith         s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
9558865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n20] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */
9568865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
95747c6ae99SBarry Smith       }
95847c6ae99SBarry Smith     }
95947c6ae99SBarry Smith 
96047c6ae99SBarry Smith     for (i=0; i<y; i++) {
96147c6ae99SBarry Smith       if (n21 >= 0) { /* directly left */
962ce00eea3SSatish Balay         x_t = lx[n21 % m];
96347c6ae99SBarry Smith         y_t = y;
96447c6ae99SBarry Smith         /* z_t = lz[n21 / (m*n)]; */
96547c6ae99SBarry Smith         s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
9668865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n21] + (i+1)*x_t - s_x;  /* 2d case */
9678865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
96847c6ae99SBarry Smith       }
96947c6ae99SBarry Smith 
97047c6ae99SBarry Smith       if (n22 >= 0) { /* middle */
97147c6ae99SBarry Smith         x_t = x;
97247c6ae99SBarry Smith         y_t = y;
97347c6ae99SBarry Smith         /* z_t = lz[n22 / (m*n)]; */
97447c6ae99SBarry Smith         s_t = bases[n22] + i*x_t + k*x_t*y_t;
9758865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n22] + i*x_t; /* 2d case */
9768865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
977bff4a2f0SMatthew G. Knepley       } else if (bz == DM_BOUNDARY_MIRROR) {
9788865f1eaSKarl Rupp         for (j=0; j<x; j++) idx[nn++] = 0;
97947c6ae99SBarry Smith       }
98047c6ae99SBarry Smith 
98147c6ae99SBarry Smith       if (n23 >= 0) { /* directly right */
982ce00eea3SSatish Balay         x_t = lx[n23 % m];
98347c6ae99SBarry Smith         y_t = y;
98447c6ae99SBarry Smith         /* z_t = lz[n23 / (m*n)]; */
98547c6ae99SBarry Smith         s_t = bases[n23] + i*x_t + k*x_t*y_t;
9868865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n23] + i*x_t; /* 2d case */
9878865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
98847c6ae99SBarry Smith       }
98947c6ae99SBarry Smith     }
99047c6ae99SBarry Smith 
99147c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
99247c6ae99SBarry Smith       if (n24 >= 0) { /* left above */
993ce00eea3SSatish Balay         x_t = lx[n24 % m];
99447c6ae99SBarry Smith         y_t = ly[(n24 % (m*n))/m];
99547c6ae99SBarry Smith         /* z_t = lz[n24 / (m*n)]; */
99647c6ae99SBarry Smith         s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
9978865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n24] + i*x_t - s_x; /* 2d case */
9988865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
99947c6ae99SBarry Smith       }
100047c6ae99SBarry Smith       if (n25 >= 0) { /* directly above */
100147c6ae99SBarry Smith         x_t = x;
100247c6ae99SBarry Smith         y_t = ly[(n25 % (m*n))/m];
100347c6ae99SBarry Smith         /* z_t = lz[n25 / (m*n)]; */
100447c6ae99SBarry Smith         s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
10058865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n25] + (i-1)*x_t; /* 2d case */
10068865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
100747c6ae99SBarry Smith       }
100847c6ae99SBarry Smith       if (n26 >= 0) { /* right above */
1009ce00eea3SSatish Balay         x_t = lx[n26 % m];
101047c6ae99SBarry Smith         y_t = ly[(n26 % (m*n))/m];
101147c6ae99SBarry Smith         /* z_t = lz[n26 / (m*n)]; */
101247c6ae99SBarry Smith         s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
10138865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n26] + (i-1)*x_t; /* 2d case */
10148865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
101547c6ae99SBarry Smith       }
101647c6ae99SBarry Smith     }
101747c6ae99SBarry Smith   }
1018ce00eea3SSatish Balay 
1019b1fb7eb7SBarry Smith   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from);CHKERRQ(ierr);
102047c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
10213bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr);
1022fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
1023fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
102447c6ae99SBarry Smith 
1025d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
102647c6ae99SBarry Smith     n0  = sn0;  n1  = sn1;  n2  = sn2;  n3  = sn3;  n5  = sn5;  n6  = sn6; n7 = sn7;
102747c6ae99SBarry Smith     n8  = sn8;  n9  = sn9;  n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18;
102847c6ae99SBarry Smith     n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25;
102947c6ae99SBarry Smith     n26 = sn26;
1030ce00eea3SSatish Balay   }
103147c6ae99SBarry Smith 
1032288e7d53SBarry Smith   if (((stencil_type == DMDA_STENCIL_STAR) || (bx != DM_BOUNDARY_PERIODIC && bx) || (by != DM_BOUNDARY_PERIODIC && by) || (bz != DM_BOUNDARY_PERIODIC && bz))) {
1033ce00eea3SSatish Balay     /*
1034ce00eea3SSatish Balay         Recompute the local to global mappings, this time keeping the
1035ce00eea3SSatish Balay       information about the cross corner processor numbers.
1036ce00eea3SSatish Balay     */
103747c6ae99SBarry Smith     nn = 0;
103847c6ae99SBarry Smith     /* Bottom Level */
103947c6ae99SBarry Smith     for (k=0; k<s_z; k++) {
104047c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
104147c6ae99SBarry Smith         if (n0 >= 0) { /* left below */
1042ce00eea3SSatish Balay           x_t = lx[n0 % m];
104347c6ae99SBarry Smith           y_t = ly[(n0 % (m*n))/m];
104447c6ae99SBarry Smith           z_t = lz[n0 / (m*n)];
104547c6ae99SBarry 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;
10468865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1047ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Ys-ys < 0 && Zs-zs < 0) {
10488865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
104947c6ae99SBarry Smith         }
105047c6ae99SBarry Smith         if (n1 >= 0) { /* directly below */
105147c6ae99SBarry Smith           x_t = x;
105247c6ae99SBarry Smith           y_t = ly[(n1 % (m*n))/m];
105347c6ae99SBarry Smith           z_t = lz[n1 / (m*n)];
105447c6ae99SBarry 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;
10558865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1056ce00eea3SSatish Balay         } else if (Ys-ys < 0 && Zs-zs < 0) {
10578865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
105847c6ae99SBarry Smith         }
105947c6ae99SBarry Smith         if (n2 >= 0) { /* right below */
1060ce00eea3SSatish Balay           x_t = lx[n2 % m];
106147c6ae99SBarry Smith           y_t = ly[(n2 % (m*n))/m];
106247c6ae99SBarry Smith           z_t = lz[n2 / (m*n)];
106347c6ae99SBarry 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;
10648865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1065ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Ys-ys < 0 && Zs-zs < 0) {
10668865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
106747c6ae99SBarry Smith         }
106847c6ae99SBarry Smith       }
106947c6ae99SBarry Smith 
107047c6ae99SBarry Smith       for (i=0; i<y; i++) {
107147c6ae99SBarry Smith         if (n3 >= 0) { /* directly left */
1072ce00eea3SSatish Balay           x_t = lx[n3 % m];
107347c6ae99SBarry Smith           y_t = y;
107447c6ae99SBarry Smith           z_t = lz[n3 / (m*n)];
107547c6ae99SBarry Smith           s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
10768865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1077ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Zs-zs < 0) {
10788865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
107947c6ae99SBarry Smith         }
108047c6ae99SBarry Smith 
108147c6ae99SBarry Smith         if (n4 >= 0) { /* middle */
108247c6ae99SBarry Smith           x_t = x;
108347c6ae99SBarry Smith           y_t = y;
108447c6ae99SBarry Smith           z_t = lz[n4 / (m*n)];
108547c6ae99SBarry Smith           s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
10868865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1087ce00eea3SSatish Balay         } else if (Zs-zs < 0) {
1088bff4a2f0SMatthew G. Knepley           if (bz == DM_BOUNDARY_MIRROR) {
10898865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = 0;
1090c2859e5eSBarry Smith           } else {
10918865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = -1;
109247c6ae99SBarry Smith           }
1093c2859e5eSBarry Smith         }
109447c6ae99SBarry Smith 
109547c6ae99SBarry Smith         if (n5 >= 0) { /* directly right */
1096ce00eea3SSatish Balay           x_t = lx[n5 % m];
109747c6ae99SBarry Smith           y_t = y;
109847c6ae99SBarry Smith           z_t = lz[n5 / (m*n)];
109947c6ae99SBarry Smith           s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
11008865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1101ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Zs-zs < 0) {
11028865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
110347c6ae99SBarry Smith         }
110447c6ae99SBarry Smith       }
110547c6ae99SBarry Smith 
110647c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
110747c6ae99SBarry Smith         if (n6 >= 0) { /* left above */
1108ce00eea3SSatish Balay           x_t = lx[n6 % m];
110947c6ae99SBarry Smith           y_t = ly[(n6 % (m*n))/m];
111047c6ae99SBarry Smith           z_t = lz[n6 / (m*n)];
111147c6ae99SBarry Smith           s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
11128865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1113ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ye-Ye < 0 && Zs-zs < 0) {
11148865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
111547c6ae99SBarry Smith         }
111647c6ae99SBarry Smith         if (n7 >= 0) { /* directly above */
111747c6ae99SBarry Smith           x_t = x;
111847c6ae99SBarry Smith           y_t = ly[(n7 % (m*n))/m];
111947c6ae99SBarry Smith           z_t = lz[n7 / (m*n)];
112047c6ae99SBarry Smith           s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
11218865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1122ce00eea3SSatish Balay         } else if (ye-Ye < 0 && Zs-zs < 0) {
11238865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
112447c6ae99SBarry Smith         }
112547c6ae99SBarry Smith         if (n8 >= 0) { /* right above */
1126ce00eea3SSatish Balay           x_t = lx[n8 % m];
112747c6ae99SBarry Smith           y_t = ly[(n8 % (m*n))/m];
112847c6ae99SBarry Smith           z_t = lz[n8 / (m*n)];
112947c6ae99SBarry Smith           s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
11308865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1131ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ye-Ye < 0 && Zs-zs < 0) {
11328865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
113347c6ae99SBarry Smith         }
113447c6ae99SBarry Smith       }
113547c6ae99SBarry Smith     }
113647c6ae99SBarry Smith 
113747c6ae99SBarry Smith     /* Middle Level */
113847c6ae99SBarry Smith     for (k=0; k<z; k++) {
113947c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
114047c6ae99SBarry Smith         if (n9 >= 0) { /* left below */
1141ce00eea3SSatish Balay           x_t = lx[n9 % m];
114247c6ae99SBarry Smith           y_t = ly[(n9 % (m*n))/m];
114347c6ae99SBarry Smith           /* z_t = z; */
114447c6ae99SBarry Smith           s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
11458865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1146ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Ys-ys < 0) {
11478865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
114847c6ae99SBarry Smith         }
114947c6ae99SBarry Smith         if (n10 >= 0) { /* directly below */
115047c6ae99SBarry Smith           x_t = x;
115147c6ae99SBarry Smith           y_t = ly[(n10 % (m*n))/m];
115247c6ae99SBarry Smith           /* z_t = z; */
115347c6ae99SBarry Smith           s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
11548865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1155ce00eea3SSatish Balay         } else if (Ys-ys < 0) {
1156bff4a2f0SMatthew G. Knepley           if (by == DM_BOUNDARY_MIRROR) {
11578865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = -1;
1158c2859e5eSBarry Smith           } else {
11598865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = -1;
1160c2859e5eSBarry Smith           }
116147c6ae99SBarry Smith         }
116247c6ae99SBarry Smith         if (n11 >= 0) { /* right below */
1163ce00eea3SSatish Balay           x_t = lx[n11 % m];
116447c6ae99SBarry Smith           y_t = ly[(n11 % (m*n))/m];
116547c6ae99SBarry Smith           /* z_t = z; */
116647c6ae99SBarry Smith           s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
11678865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1168ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Ys-ys < 0) {
11698865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
117047c6ae99SBarry Smith         }
117147c6ae99SBarry Smith       }
117247c6ae99SBarry Smith 
117347c6ae99SBarry Smith       for (i=0; i<y; i++) {
117447c6ae99SBarry Smith         if (n12 >= 0) { /* directly left */
1175ce00eea3SSatish Balay           x_t = lx[n12 % m];
117647c6ae99SBarry Smith           y_t = y;
117747c6ae99SBarry Smith           /* z_t = z; */
117847c6ae99SBarry Smith           s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
11798865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1180ce00eea3SSatish Balay         } else if (Xs-xs < 0) {
1181bff4a2f0SMatthew G. Knepley           if (bx == DM_BOUNDARY_MIRROR) {
11828865f1eaSKarl Rupp             for (j=0; j<s_x; j++) idx[nn++] = 0;
1183c2859e5eSBarry Smith           } else {
11848865f1eaSKarl Rupp             for (j=0; j<s_x; j++) idx[nn++] = -1;
118547c6ae99SBarry Smith           }
1186c2859e5eSBarry Smith         }
118747c6ae99SBarry Smith 
118847c6ae99SBarry Smith         /* Interior */
118947c6ae99SBarry Smith         s_t = bases[rank] + i*x + k*x*y;
11908865f1eaSKarl Rupp         for (j=0; j<x; j++) idx[nn++] = s_t++;
119147c6ae99SBarry Smith 
119247c6ae99SBarry Smith         if (n14 >= 0) { /* directly right */
1193ce00eea3SSatish Balay           x_t = lx[n14 % m];
119447c6ae99SBarry Smith           y_t = y;
119547c6ae99SBarry Smith           /* z_t = z; */
119647c6ae99SBarry Smith           s_t = bases[n14] + i*x_t + k*x_t*y_t;
11978865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1198ce00eea3SSatish Balay         } else if (xe-Xe < 0) {
1199bff4a2f0SMatthew G. Knepley           if (bx == DM_BOUNDARY_MIRROR) {
12008865f1eaSKarl Rupp             for (j=0; j<s_x; j++) idx[nn++] = 0;
1201c2859e5eSBarry Smith           } else {
12028865f1eaSKarl Rupp             for (j=0; j<s_x; j++) idx[nn++] = -1;
120347c6ae99SBarry Smith           }
120447c6ae99SBarry Smith         }
1205c2859e5eSBarry Smith       }
120647c6ae99SBarry Smith 
120747c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
120847c6ae99SBarry Smith         if (n15 >= 0) { /* left above */
1209ce00eea3SSatish Balay           x_t = lx[n15 % m];
121047c6ae99SBarry Smith           y_t = ly[(n15 % (m*n))/m];
121147c6ae99SBarry Smith           /* z_t = z; */
121247c6ae99SBarry Smith           s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
12138865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1214ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ye-Ye < 0) {
12158865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
121647c6ae99SBarry Smith         }
121747c6ae99SBarry Smith         if (n16 >= 0) { /* directly above */
121847c6ae99SBarry Smith           x_t = x;
121947c6ae99SBarry Smith           y_t = ly[(n16 % (m*n))/m];
122047c6ae99SBarry Smith           /* z_t = z; */
122147c6ae99SBarry Smith           s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
12228865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1223ce00eea3SSatish Balay         } else if (ye-Ye < 0) {
1224bff4a2f0SMatthew G. Knepley           if (by == DM_BOUNDARY_MIRROR) {
12258865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = 0;
1226c2859e5eSBarry Smith           } else {
12278865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = -1;
122847c6ae99SBarry Smith           }
1229c2859e5eSBarry Smith         }
123047c6ae99SBarry Smith         if (n17 >= 0) { /* right above */
1231ce00eea3SSatish Balay           x_t = lx[n17 % m];
123247c6ae99SBarry Smith           y_t = ly[(n17 % (m*n))/m];
123347c6ae99SBarry Smith           /* z_t = z; */
123447c6ae99SBarry Smith           s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
12358865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1236ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ye-Ye < 0) {
12378865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
123847c6ae99SBarry Smith         }
123947c6ae99SBarry Smith       }
124047c6ae99SBarry Smith     }
124147c6ae99SBarry Smith 
124247c6ae99SBarry Smith     /* Upper Level */
124347c6ae99SBarry Smith     for (k=0; k<s_z; k++) {
124447c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
124547c6ae99SBarry Smith         if (n18 >= 0) { /* left below */
1246ce00eea3SSatish Balay           x_t = lx[n18 % m];
124747c6ae99SBarry Smith           y_t = ly[(n18 % (m*n))/m];
124847c6ae99SBarry Smith           /* z_t = lz[n18 / (m*n)]; */
124947c6ae99SBarry Smith           s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
12508865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1251ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Ys-ys < 0 && ze-Ze < 0) {
12528865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
125347c6ae99SBarry Smith         }
125447c6ae99SBarry Smith         if (n19 >= 0) { /* directly below */
125547c6ae99SBarry Smith           x_t = x;
125647c6ae99SBarry Smith           y_t = ly[(n19 % (m*n))/m];
125747c6ae99SBarry Smith           /* z_t = lz[n19 / (m*n)]; */
125847c6ae99SBarry Smith           s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
12598865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1260ce00eea3SSatish Balay         } else if (Ys-ys < 0 && ze-Ze < 0) {
12618865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
126247c6ae99SBarry Smith         }
126347c6ae99SBarry Smith         if (n20 >= 0) { /* right below */
1264ce00eea3SSatish Balay           x_t = lx[n20 % m];
126547c6ae99SBarry Smith           y_t = ly[(n20 % (m*n))/m];
126647c6ae99SBarry Smith           /* z_t = lz[n20 / (m*n)]; */
126747c6ae99SBarry Smith           s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
12688865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1269ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Ys-ys < 0 && ze-Ze < 0) {
12708865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
127147c6ae99SBarry Smith         }
127247c6ae99SBarry Smith       }
127347c6ae99SBarry Smith 
127447c6ae99SBarry Smith       for (i=0; i<y; i++) {
127547c6ae99SBarry Smith         if (n21 >= 0) { /* directly left */
1276ce00eea3SSatish Balay           x_t = lx[n21 % m];
127747c6ae99SBarry Smith           y_t = y;
127847c6ae99SBarry Smith           /* z_t = lz[n21 / (m*n)]; */
127947c6ae99SBarry Smith           s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
12808865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1281ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ze-Ze < 0) {
12828865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
128347c6ae99SBarry Smith         }
128447c6ae99SBarry Smith 
128547c6ae99SBarry Smith         if (n22 >= 0) { /* middle */
128647c6ae99SBarry Smith           x_t = x;
128747c6ae99SBarry Smith           y_t = y;
128847c6ae99SBarry Smith           /* z_t = lz[n22 / (m*n)]; */
128947c6ae99SBarry Smith           s_t = bases[n22] + i*x_t + k*x_t*y_t;
12908865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1291ce00eea3SSatish Balay         } else if (ze-Ze < 0) {
1292bff4a2f0SMatthew G. Knepley           if (bz == DM_BOUNDARY_MIRROR) {
12938865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = 0;
1294c2859e5eSBarry Smith           } else {
12958865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = -1;
129647c6ae99SBarry Smith           }
1297c2859e5eSBarry Smith         }
129847c6ae99SBarry Smith 
129947c6ae99SBarry Smith         if (n23 >= 0) { /* directly right */
1300ce00eea3SSatish Balay           x_t = lx[n23 % m];
130147c6ae99SBarry Smith           y_t = y;
130247c6ae99SBarry Smith           /* z_t = lz[n23 / (m*n)]; */
130347c6ae99SBarry Smith           s_t = bases[n23] + i*x_t + k*x_t*y_t;
13048865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1305ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ze-Ze < 0) {
13068865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
130747c6ae99SBarry Smith         }
130847c6ae99SBarry Smith       }
130947c6ae99SBarry Smith 
131047c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
131147c6ae99SBarry Smith         if (n24 >= 0) { /* left above */
1312ce00eea3SSatish Balay           x_t = lx[n24 % m];
131347c6ae99SBarry Smith           y_t = ly[(n24 % (m*n))/m];
131447c6ae99SBarry Smith           /* z_t = lz[n24 / (m*n)]; */
131547c6ae99SBarry Smith           s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
13168865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1317ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ye-Ye < 0 && ze-Ze < 0) {
13188865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
131947c6ae99SBarry Smith         }
132047c6ae99SBarry Smith         if (n25 >= 0) { /* directly above */
132147c6ae99SBarry Smith           x_t = x;
132247c6ae99SBarry Smith           y_t = ly[(n25 % (m*n))/m];
132347c6ae99SBarry Smith           /* z_t = lz[n25 / (m*n)]; */
132447c6ae99SBarry Smith           s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
13258865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1326ce00eea3SSatish Balay         } else if (ye-Ye < 0 && ze-Ze < 0) {
13278865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
132847c6ae99SBarry Smith         }
132947c6ae99SBarry Smith         if (n26 >= 0) { /* right above */
1330ce00eea3SSatish Balay           x_t = lx[n26 % m];
133147c6ae99SBarry Smith           y_t = ly[(n26 % (m*n))/m];
133247c6ae99SBarry Smith           /* z_t = lz[n26 / (m*n)]; */
133347c6ae99SBarry Smith           s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
13348865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1335ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ye-Ye < 0 && ze-Ze < 0) {
13368865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
133747c6ae99SBarry Smith         }
133847c6ae99SBarry Smith       }
133947c6ae99SBarry Smith     }
134047c6ae99SBarry Smith   }
134147c6ae99SBarry Smith   /*
134247c6ae99SBarry Smith      Set the local to global ordering in the global vector, this allows use
134347c6ae99SBarry Smith      of VecSetValuesLocal().
134447c6ae99SBarry Smith   */
134545b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr);
13463bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr);
134747c6ae99SBarry Smith 
1348ce00eea3SSatish Balay   ierr  = PetscFree2(bases,ldims);CHKERRQ(ierr);
1349ce00eea3SSatish Balay   dd->m = m;  dd->n  = n;  dd->p  = p;
1350ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
1351ce00eea3SSatish Balay   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = zs; dd->ze = ze;
1352ce00eea3SSatish Balay   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = Zs; dd->Ze = Ze;
1353ce00eea3SSatish Balay 
1354fcfd50ebSBarry Smith   ierr = VecDestroy(&local);CHKERRQ(ierr);
1355fcfd50ebSBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
1356ce00eea3SSatish Balay 
1357ce00eea3SSatish Balay   dd->gtol      = gtol;
1358ce00eea3SSatish Balay   dd->base      = base;
1359ce00eea3SSatish Balay   da->ops->view = DMView_DA_3d;
13600298fd71SBarry Smith   dd->ltol      = NULL;
13610298fd71SBarry Smith   dd->ao        = NULL;
136247c6ae99SBarry Smith   PetscFunctionReturn(0);
136347c6ae99SBarry Smith }
1364564755cdSBarry Smith 
136547c6ae99SBarry Smith 
136647c6ae99SBarry Smith /*@C
1367aa219208SBarry Smith    DMDACreate3d - Creates an object that will manage the communication of three-dimensional
136847c6ae99SBarry Smith    regular array data that is distributed across some processors.
136947c6ae99SBarry Smith 
137047c6ae99SBarry Smith    Collective on MPI_Comm
137147c6ae99SBarry Smith 
137247c6ae99SBarry Smith    Input Parameters:
137347c6ae99SBarry Smith +  comm - MPI communicator
13741321219cSEthan Coon .  bx,by,bz - type of ghost nodes the array have.
1375bff4a2f0SMatthew G. Knepley          Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC.
1376aa219208SBarry Smith .  stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX)
1377897f7067SBarry Smith .  M,N,P - global dimension in each direction of the array
137847c6ae99SBarry Smith .  m,n,p - corresponding number of processors in each dimension
137947c6ae99SBarry Smith            (or PETSC_DECIDE to have calculated)
138047c6ae99SBarry Smith .  dof - number of degrees of freedom per node
138110d7c030SMatthew G Knepley .  s - stencil width
138210d7c030SMatthew G Knepley -  lx, ly, lz - arrays containing the number of nodes in each cell along
13830298fd71SBarry Smith           the x, y, and z coordinates, or NULL. If non-null, these
138447c6ae99SBarry Smith           must be of length as m,n,p and the corresponding
138547c6ae99SBarry Smith           m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of
138647c6ae99SBarry Smith           the ly[] must N, sum of the lz[] must be P
138747c6ae99SBarry Smith 
138847c6ae99SBarry Smith    Output Parameter:
138947c6ae99SBarry Smith .  da - the resulting distributed array object
139047c6ae99SBarry Smith 
139147c6ae99SBarry Smith    Options Database Key:
1392706a11cbSBarry Smith +  -dm_view - Calls DMView() at the conclusion of DMDACreate3d()
139347c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
139447c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
139547c6ae99SBarry Smith .  -da_grid_z <nz> - number of grid points in z direction, if P < 0
139647c6ae99SBarry Smith .  -da_processors_x <MX> - number of processors in x direction
139747c6ae99SBarry Smith .  -da_processors_y <MY> - number of processors in y direction
139847c6ae99SBarry Smith .  -da_processors_z <MZ> - number of processors in z direction
1399e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement ratio in x direction
1400e0f5d30fSBarry Smith .  -da_refine_y <ry> - refinement ratio in y direction
1401e0f5d30fSBarry Smith .  -da_refine_z <rz>- refinement ratio in z directio
1402e0f5d30fSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating it, , if M, N, or P < 0
140347c6ae99SBarry Smith 
140447c6ae99SBarry Smith    Level: beginner
140547c6ae99SBarry Smith 
140647c6ae99SBarry Smith    Notes:
1407aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
1408aa219208SBarry Smith    standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
140947c6ae99SBarry Smith    the standard 27-pt stencil.
141047c6ae99SBarry Smith 
1411aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
1412564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
1413564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
141447c6ae99SBarry Smith 
1415897f7067SBarry Smith    You must call DMSetUp() after this call before using this DM.
1416897f7067SBarry Smith 
1417897f7067SBarry Smith    If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call
1418897f7067SBarry Smith    but before DMSetUp().
1419897f7067SBarry Smith 
142047c6ae99SBarry Smith .keywords: distributed array, create, three-dimensional
142147c6ae99SBarry Smith 
1422aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate2d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
142399f0b487SRichard Tran Mills           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(),
1424d461ba97SBarry Smith           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
142547c6ae99SBarry Smith 
142647c6ae99SBarry Smith @*/
1427bff4a2f0SMatthew G. Knepley PetscErrorCode  DMDACreate3d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMBoundaryType bz,DMDAStencilType stencil_type,PetscInt M,
14289a42bb27SBarry 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)
142947c6ae99SBarry Smith {
143047c6ae99SBarry Smith   PetscErrorCode ierr;
143147c6ae99SBarry Smith 
143247c6ae99SBarry Smith   PetscFunctionBegin;
1433aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
1434c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*da, 3);CHKERRQ(ierr);
1435aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, P);CHKERRQ(ierr);
1436aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, p);CHKERRQ(ierr);
14371321219cSEthan Coon   ierr = DMDASetBoundaryType(*da, bx, by, bz);CHKERRQ(ierr);
1438aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
1439aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
1440aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
1441aa219208SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, lz);CHKERRQ(ierr);
144247c6ae99SBarry Smith   PetscFunctionReturn(0);
144347c6ae99SBarry Smith }
1444