xref: /petsc/src/dm/impls/da/da3.c (revision 6f951b9505d1f14bccad7423b8c71a09c758204a)
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 
7c6db04a5SJed Brown #include <private/daimpl.h>     /*I   "petscdmda.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 
337b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
3447c6ae99SBarry Smith     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
3547c6ae99SBarry Smith     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) {
36aa219208SBarry Smith       DMDALocalInfo info;
37aa219208SBarry Smith       ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
3847c6ae99SBarry 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);
3947c6ae99SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D, Y range of indices: %D %D, Z range of indices: %D %D\n",
4047c6ae99SBarry Smith                                                 info.xs,info.xs+info.xm,info.ys,info.ys+info.ym,info.zs,info.zs+info.zm);CHKERRQ(ierr);
4147c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX)
4247c6ae99SBarry Smith       if (dd->coordinates) {
4347c6ae99SBarry Smith         PetscInt        last;
4447c6ae99SBarry Smith         const PetscReal *coors;
4547c6ae99SBarry Smith         ierr = VecGetArrayRead(dd->coordinates,&coors);CHKERRQ(ierr);
4647c6ae99SBarry Smith         ierr = VecGetLocalSize(dd->coordinates,&last);CHKERRQ(ierr);
4747c6ae99SBarry Smith         last = last - 3;
4847c6ae99SBarry 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);
4947c6ae99SBarry Smith         ierr = VecRestoreArrayRead(dd->coordinates,&coors);CHKERRQ(ierr);
5047c6ae99SBarry Smith       }
5147c6ae99SBarry Smith #endif
5247c6ae99SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
537b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
54616157d6SJed Brown     } else {
55616157d6SJed Brown       ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr);
5647c6ae99SBarry Smith     }
5747c6ae99SBarry Smith   } else if (isdraw) {
5847c6ae99SBarry Smith     PetscDraw       draw;
5947c6ae99SBarry Smith     PetscReal     ymin = -1.0,ymax = (PetscReal)dd->N;
6047c6ae99SBarry Smith     PetscReal     xmin = -1.0,xmax = (PetscReal)((dd->M+2)*dd->P),x,y,ycoord,xcoord;
6147c6ae99SBarry Smith     PetscInt        k,plane,base,*idx;
6247c6ae99SBarry Smith     char       node[10];
6347c6ae99SBarry Smith     PetscBool  isnull;
6447c6ae99SBarry Smith 
6547c6ae99SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
6647c6ae99SBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
6747c6ae99SBarry Smith     ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
6847c6ae99SBarry Smith     ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
6947c6ae99SBarry Smith 
7047c6ae99SBarry Smith     /* first processor draw all node lines */
7147c6ae99SBarry Smith     if (!rank) {
7247c6ae99SBarry Smith       for (k=0; k<dd->P; k++) {
7347c6ae99SBarry Smith         ymin = 0.0; ymax = (PetscReal)(dd->N - 1);
7447c6ae99SBarry Smith         for (xmin=(PetscReal)(k*(dd->M+1)); xmin<(PetscReal)(dd->M+(k*(dd->M+1))); xmin++) {
7547c6ae99SBarry Smith           ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
7647c6ae99SBarry Smith         }
7747c6ae99SBarry Smith 
7847c6ae99SBarry Smith         xmin = (PetscReal)(k*(dd->M+1)); xmax = xmin + (PetscReal)(dd->M - 1);
7947c6ae99SBarry Smith         for (ymin=0; ymin<(PetscReal)dd->N; ymin++) {
8047c6ae99SBarry Smith           ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
8147c6ae99SBarry Smith         }
8247c6ae99SBarry Smith       }
8347c6ae99SBarry Smith     }
8447c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
8547c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
8647c6ae99SBarry Smith 
8747c6ae99SBarry Smith     for (k=0; k<dd->P; k++) {  /*Go through and draw for each plane*/
8847c6ae99SBarry Smith       if ((k >= dd->zs) && (k < dd->ze)) {
8947c6ae99SBarry Smith         /* draw my box */
9047c6ae99SBarry Smith         ymin = dd->ys;
9147c6ae99SBarry Smith         ymax = dd->ye - 1;
9247c6ae99SBarry Smith         xmin = dd->xs/dd->w    + (dd->M+1)*k;
9347c6ae99SBarry Smith         xmax =(dd->xe-1)/dd->w + (dd->M+1)*k;
9447c6ae99SBarry Smith 
9547c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
9647c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
9747c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
9847c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
9947c6ae99SBarry Smith 
10047c6ae99SBarry Smith         xmin = dd->xs/dd->w;
10147c6ae99SBarry Smith         xmax =(dd->xe-1)/dd->w;
10247c6ae99SBarry Smith 
10347c6ae99SBarry Smith         /* put in numbers*/
10447c6ae99SBarry Smith         base = (dd->base+(dd->xe-dd->xs)*(dd->ye-dd->ys)*(k-dd->zs))/dd->w;
10547c6ae99SBarry Smith 
10647c6ae99SBarry Smith         /* Identify which processor owns the box */
10747c6ae99SBarry Smith         sprintf(node,"%d",rank);
10847c6ae99SBarry Smith         ierr = PetscDrawString(draw,xmin+(dd->M+1)*k+.2,ymin+.3,PETSC_DRAW_RED,node);CHKERRQ(ierr);
10947c6ae99SBarry Smith 
11047c6ae99SBarry Smith         for (y=ymin; y<=ymax; y++) {
11147c6ae99SBarry Smith           for (x=xmin+(dd->M+1)*k; x<=xmax+(dd->M+1)*k; x++) {
11247c6ae99SBarry Smith             sprintf(node,"%d",(int)base++);
11347c6ae99SBarry Smith             ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
11447c6ae99SBarry Smith           }
11547c6ae99SBarry Smith         }
11647c6ae99SBarry Smith 
11747c6ae99SBarry Smith       }
11847c6ae99SBarry Smith     }
11947c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
12047c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
12147c6ae99SBarry Smith 
12247c6ae99SBarry Smith     for (k=0-dd->s; k<dd->P+dd->s; k++) {
12347c6ae99SBarry Smith       /* Go through and draw for each plane */
12447c6ae99SBarry Smith       if ((k >= dd->Zs) && (k < dd->Ze)) {
12547c6ae99SBarry Smith 
12647c6ae99SBarry Smith         /* overlay ghost numbers, useful for error checking */
12747c6ae99SBarry Smith         base = (dd->Xe-dd->Xs)*(dd->Ye-dd->Ys)*(k-dd->Zs); idx = dd->idx;
12847c6ae99SBarry Smith         plane=k;
12947c6ae99SBarry Smith         /* Keep z wrap around points on the dradrawg */
13047c6ae99SBarry Smith         if (k<0)    { plane=dd->P+k; }
13147c6ae99SBarry Smith         if (k>=dd->P) { plane=k-dd->P; }
13247c6ae99SBarry Smith         ymin = dd->Ys; ymax = dd->Ye;
13347c6ae99SBarry Smith         xmin = (dd->M+1)*plane*dd->w;
13447c6ae99SBarry Smith         xmax = (dd->M+1)*plane*dd->w+dd->M*dd->w;
13547c6ae99SBarry Smith         for (y=ymin; y<ymax; y++) {
13647c6ae99SBarry Smith           for (x=xmin+dd->Xs; x<xmin+dd->Xe; x+=dd->w) {
13747c6ae99SBarry Smith             sprintf(node,"%d",(int)(idx[base]/dd->w));
13847c6ae99SBarry Smith             ycoord = y;
13947c6ae99SBarry Smith             /*Keep y wrap around points on drawing */
14047c6ae99SBarry Smith             if (y<0)      { ycoord = dd->N+y; }
14147c6ae99SBarry Smith 
14247c6ae99SBarry Smith             if (y>=dd->N) { ycoord = y-dd->N; }
14347c6ae99SBarry Smith             xcoord = x;   /* Keep x wrap points on drawing */
14447c6ae99SBarry Smith 
14547c6ae99SBarry Smith             if (x<xmin)  { xcoord = xmax - (xmin-x); }
14647c6ae99SBarry Smith             if (x>=xmax) { xcoord = xmin + (x-xmax); }
14747c6ae99SBarry Smith             ierr = PetscDrawString(draw,xcoord/dd->w,ycoord,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
14847c6ae99SBarry Smith             base+=dd->w;
14947c6ae99SBarry Smith           }
15047c6ae99SBarry Smith         }
15147c6ae99SBarry Smith       }
15247c6ae99SBarry Smith     }
15347c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
15447c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
1559a42bb27SBarry Smith   } else if (isbinary){
1569a42bb27SBarry Smith     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
1579a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1589a42bb27SBarry Smith   } else if (ismatlab) {
1599a42bb27SBarry Smith     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
1609a42bb27SBarry Smith #endif
161aa219208SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for DMDA 1d",((PetscObject)viewer)->type_name);
16247c6ae99SBarry Smith   PetscFunctionReturn(0);
16347c6ae99SBarry Smith }
16447c6ae99SBarry Smith 
16547c6ae99SBarry Smith #undef __FUNCT__
1669a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_3D"
1677087cfbeSBarry Smith PetscErrorCode  DMSetUp_DA_3D(DM da)
16847c6ae99SBarry Smith {
16947c6ae99SBarry Smith   DM_DA                  *dd           = (DM_DA*)da->data;
17047c6ae99SBarry Smith   const PetscInt         M            = dd->M;
17147c6ae99SBarry Smith   const PetscInt         N            = dd->N;
17247c6ae99SBarry Smith   const PetscInt         P            = dd->P;
17347c6ae99SBarry Smith   PetscInt               m            = dd->m;
17447c6ae99SBarry Smith   PetscInt               n            = dd->n;
17547c6ae99SBarry Smith   PetscInt               p            = dd->p;
17647c6ae99SBarry Smith   const PetscInt         dof          = dd->w;
17747c6ae99SBarry Smith   const PetscInt         s            = dd->s;
1781321219cSEthan Coon   const DMDABoundaryType bx         = dd->bx;
1791321219cSEthan Coon   const DMDABoundaryType by         = dd->by;
1801321219cSEthan Coon   const DMDABoundaryType bz         = dd->bz;
181aa219208SBarry Smith   const DMDAStencilType  stencil_type = dd->stencil_type;
18247c6ae99SBarry Smith   PetscInt               *lx           = dd->lx;
18347c6ae99SBarry Smith   PetscInt               *ly           = dd->ly;
18447c6ae99SBarry Smith   PetscInt               *lz           = dd->lz;
18547c6ae99SBarry Smith   MPI_Comm               comm;
18647c6ae99SBarry Smith   PetscMPIInt            rank,size;
187ce00eea3SSatish Balay   PetscInt               xs = 0,xe,ys = 0,ye,zs = 0,ze,x = 0,y = 0,z = 0;
188ce00eea3SSatish Balay   PetscInt               Xs,Xe,Ys,Ye,Zs,Ze,IXs,IXe,IYs,IYe,IZs,IZe,start,end,pm;
189ce00eea3SSatish Balay   PetscInt               left,right,up,down,bottom,top,i,j,k,*idx,*idx_cpy,nn;
190ce00eea3SSatish Balay   const PetscInt         *idx_full;
19147c6ae99SBarry Smith   PetscInt               n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14;
19247c6ae99SBarry Smith   PetscInt               n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26;
193db87c5ecSEthan Coon   PetscInt               *bases,*ldims,base,x_t,y_t,z_t,s_t,count,s_x,s_y,s_z;
19447c6ae99SBarry Smith   PetscInt               sn0 = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0;
19547c6ae99SBarry Smith   PetscInt               sn8 = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0;
19647c6ae99SBarry Smith   PetscInt               sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0;
19747c6ae99SBarry Smith   Vec                    local,global;
19847c6ae99SBarry Smith   VecScatter             ltog,gtol;
199ce00eea3SSatish Balay   IS                     to,from,ltogis;
200*6f951b95Secoon   PetscBool              twod;
20147c6ae99SBarry Smith   PetscErrorCode         ierr;
20247c6ae99SBarry Smith 
203*6f951b95Secoon 
20447c6ae99SBarry Smith   PetscFunctionBegin;
20547c6ae99SBarry Smith   if (dof < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
20647c6ae99SBarry Smith   if (s < 0) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
20747c6ae99SBarry Smith 
20847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr);
20947c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
21047c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
21147c6ae99SBarry Smith 
21247c6ae99SBarry Smith   dd->dim = 3;
21347c6ae99SBarry Smith   ierr = PetscMalloc(dof*sizeof(char*),&dd->fieldname);CHKERRQ(ierr);
21447c6ae99SBarry Smith   ierr = PetscMemzero(dd->fieldname,dof*sizeof(char*));CHKERRQ(ierr);
21547c6ae99SBarry Smith 
21647c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
21747c6ae99SBarry Smith     if (m < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
21847c6ae99SBarry Smith     else if (m > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
21947c6ae99SBarry Smith   }
22047c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
22147c6ae99SBarry Smith     if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
22247c6ae99SBarry Smith     else if (n > size)  SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
22347c6ae99SBarry Smith   }
22447c6ae99SBarry Smith   if (p != PETSC_DECIDE) {
22547c6ae99SBarry Smith     if (p < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %D",p);
22647c6ae99SBarry Smith     else if (p > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %D %d",p,size);
22747c6ae99SBarry Smith   }
2280716a85fSBarry 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);
22947c6ae99SBarry Smith 
23047c6ae99SBarry Smith   /* Partition the array among the processors */
23147c6ae99SBarry Smith   if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) {
23247c6ae99SBarry Smith     m = size/(n*p);
23347c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
23447c6ae99SBarry Smith     n = size/(m*p);
23547c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
23647c6ae99SBarry Smith     p = size/(m*n);
23747c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
23847c6ae99SBarry Smith     /* try for squarish distribution */
2398f1a2a5eSBarry Smith     m = (int)(0.5 + sqrt(((double)M)*((double)size)/((double)N*p)));
24047c6ae99SBarry Smith     if (!m) m = 1;
24147c6ae99SBarry Smith     while (m > 0) {
24247c6ae99SBarry Smith       n = size/(m*p);
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 p value: p = %D",p);
24747c6ae99SBarry Smith     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
24847c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
24947c6ae99SBarry Smith     /* try for squarish distribution */
2508f1a2a5eSBarry Smith     m = (int)(0.5 + sqrt(((double)M)*((double)size)/((double)P*n)));
25147c6ae99SBarry Smith     if (!m) m = 1;
25247c6ae99SBarry Smith     while (m > 0) {
25347c6ae99SBarry Smith       p = size/(m*n);
25447c6ae99SBarry Smith       if (m*n*p == size) break;
25547c6ae99SBarry Smith       m--;
25647c6ae99SBarry Smith     }
25747c6ae99SBarry Smith     if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %D",n);
25847c6ae99SBarry Smith     if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;}
25947c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
26047c6ae99SBarry Smith     /* try for squarish distribution */
2618f1a2a5eSBarry Smith     n = (int)(0.5 + sqrt(((double)N)*((double)size)/((double)P*m)));
26247c6ae99SBarry Smith     if (!n) n = 1;
26347c6ae99SBarry Smith     while (n > 0) {
26447c6ae99SBarry Smith       p = size/(m*n);
26547c6ae99SBarry Smith       if (m*n*p == size) break;
26647c6ae99SBarry Smith       n--;
26747c6ae99SBarry Smith     }
26847c6ae99SBarry Smith     if (!n) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %D",n);
26947c6ae99SBarry Smith     if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;}
27047c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
27147c6ae99SBarry Smith     /* try for squarish distribution */
2728f1a2a5eSBarry Smith     n = (PetscInt)(0.5 + pow(((double)N*N)*((double)size)/((double)P*M),(double)(1./3.)));
27347c6ae99SBarry Smith     if (!n) n = 1;
27447c6ae99SBarry Smith     while (n > 0) {
27547c6ae99SBarry Smith       pm = size/n;
27647c6ae99SBarry Smith       if (n*pm == size) break;
27747c6ae99SBarry Smith       n--;
27847c6ae99SBarry Smith     }
27947c6ae99SBarry Smith     if (!n) n = 1;
2808f1a2a5eSBarry Smith     m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)P*n)));
28147c6ae99SBarry Smith     if (!m) m = 1;
28247c6ae99SBarry Smith     while (m > 0) {
28347c6ae99SBarry Smith       p = size/(m*n);
28447c6ae99SBarry Smith       if (m*n*p == size) break;
28547c6ae99SBarry Smith       m--;
28647c6ae99SBarry Smith     }
28747c6ae99SBarry Smith     if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;}
28847c6ae99SBarry Smith   } else if (m*n*p != size) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
28947c6ae99SBarry Smith 
29047c6ae99SBarry Smith   if (m*n*p != size) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_PLIB,"Could not find good partition");
29147c6ae99SBarry Smith   if (M < m) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
29247c6ae99SBarry Smith   if (N < n) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
29347c6ae99SBarry Smith   if (P < p) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %D %D",P,p);
29447c6ae99SBarry Smith 
29547c6ae99SBarry Smith   /*
29647c6ae99SBarry Smith      Determine locally owned region
29747c6ae99SBarry Smith      [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes
29847c6ae99SBarry Smith   */
29947c6ae99SBarry Smith 
30047c6ae99SBarry Smith   if (!lx) {
30147c6ae99SBarry Smith     ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
30247c6ae99SBarry Smith     lx = dd->lx;
30347c6ae99SBarry Smith     for (i=0; i<m; i++) {
30447c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > (i % m));
30547c6ae99SBarry Smith     }
30647c6ae99SBarry Smith   }
30747c6ae99SBarry Smith   x  = lx[rank % m];
30847c6ae99SBarry Smith   xs = 0;
30947c6ae99SBarry Smith   for (i=0; i<(rank%m); i++) { xs += lx[i];}
310bcea557cSEthan Coon   if ((x < s) && ((m > 1) || (bx == DMDA_BOUNDARY_PERIODIC))) {
311bcea557cSEthan Coon     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s);
312bcea557cSEthan Coon   }
31347c6ae99SBarry Smith 
31447c6ae99SBarry Smith   if (!ly) {
31547c6ae99SBarry Smith     ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
31647c6ae99SBarry Smith     ly = dd->ly;
31747c6ae99SBarry Smith     for (i=0; i<n; i++) {
31847c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > (i % n));
31947c6ae99SBarry Smith     }
32047c6ae99SBarry Smith   }
32147c6ae99SBarry Smith   y  = ly[(rank % (m*n))/m];
322bcea557cSEthan Coon   if ((y < s) && ((n > 1) || (by == DMDA_BOUNDARY_PERIODIC))) {
323bcea557cSEthan Coon     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s);
324bcea557cSEthan Coon   }
32547c6ae99SBarry Smith   ys = 0;
32647c6ae99SBarry Smith   for (i=0; i<(rank % (m*n))/m; i++) { ys += ly[i];}
32747c6ae99SBarry Smith 
32847c6ae99SBarry Smith   if (!lz) {
32947c6ae99SBarry Smith     ierr = PetscMalloc(p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr);
33047c6ae99SBarry Smith     lz = dd->lz;
33147c6ae99SBarry Smith     for (i=0; i<p; i++) {
33247c6ae99SBarry Smith       lz[i] = P/p + ((P % p) > (i % p));
33347c6ae99SBarry Smith     }
33447c6ae99SBarry Smith   }
33547c6ae99SBarry Smith   z  = lz[rank/(m*n)];
336bcea557cSEthan Coon 
337fdc81ce8SEthan Coon   /* note this is different than x- and y-, as we will handle as an important special
338fdc81ce8SEthan Coon    case when p=P=1 and DMDA_BOUNDARY_PERIODIC and s > z.  This is to deal with 2D problems
339fdc81ce8SEthan Coon    in a 3D code.  Additional code for this case is noted with "2d case" comments */
340*6f951b95Secoon   twod = PETSC_FALSE;
341*6f951b95Secoon   if (P == 1) {
342*6f951b95Secoon     twod = PETSC_TRUE;
343*6f951b95Secoon   } else if ((z < s) && ((p > 1) || (bz == DMDA_BOUNDARY_PERIODIC))) {
344bcea557cSEthan Coon     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local z-width of domain z %D is smaller than stencil width s %D",z,s);
345bcea557cSEthan Coon   }
34647c6ae99SBarry Smith   zs = 0;
34747c6ae99SBarry Smith   for (i=0; i<(rank/(m*n)); i++) { zs += lz[i];}
34847c6ae99SBarry Smith   ye = ys + y;
34947c6ae99SBarry Smith   xe = xs + x;
35047c6ae99SBarry Smith   ze = zs + z;
35147c6ae99SBarry Smith 
35247c6ae99SBarry Smith   /* determine ghost region */
35347c6ae99SBarry Smith   /* Assume No Periodicity */
354ce00eea3SSatish Balay   if (xs-s > 0) { Xs = xs - s; IXs = xs - s; } else { Xs = 0; IXs = 0; }
355ce00eea3SSatish Balay   if (xe+s <= M) { Xe = xe + s; IXe = xe + s; } else { Xe = M; IXe = M; }
356ce00eea3SSatish Balay   if (ys-s > 0) { Ys = ys - s; IYs = ys - s; } else { Ys = 0; IYs = 0; }
357ce00eea3SSatish Balay   if (ye+s <= N) { Ye = ye + s; IYe = ye + s; } else { Ye = N; IYe = N; }
358ce00eea3SSatish Balay   if (zs-s > 0) { Zs = zs - s; IZs = zs - s; } else { Zs = 0; IZs = 0; }
359ce00eea3SSatish Balay   if (ze+s <= P) { Ze = ze + s; IZe = ze + s; } else { Ze = P; IZe = P; }
36047c6ae99SBarry Smith 
361ce00eea3SSatish Balay   /* fix for periodicity/ghosted */
3621321219cSEthan Coon   if (bx) { Xs = xs - s; Xe = xe + s; }
3631321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC) { IXs = xs - s; IXe = xe + s; }
3641321219cSEthan Coon   if (by) { Ys = ys - s; Ye = ye + s; }
3651321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC) { IYs = ys - s; IYe = ye + s; }
3661321219cSEthan Coon   if (bz) { Zs = zs - s; Ze = ze + s; }
3671321219cSEthan Coon   if (bz == DMDA_BOUNDARY_PERIODIC) { IZs = zs - s; IZe = ze + s; }
36847c6ae99SBarry Smith 
36947c6ae99SBarry Smith   /* Resize all X parameters to reflect w */
370ce00eea3SSatish Balay   s_x = s;
37147c6ae99SBarry Smith   s_y  = s;
37247c6ae99SBarry Smith   s_z  = s;
37347c6ae99SBarry Smith 
37447c6ae99SBarry Smith   /* determine starting point of each processor */
37547c6ae99SBarry Smith   nn       = x*y*z;
37647c6ae99SBarry Smith   ierr     = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr);
37747c6ae99SBarry Smith   ierr     = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
37847c6ae99SBarry Smith   bases[0] = 0;
37947c6ae99SBarry Smith   for (i=1; i<=size; i++) {
38047c6ae99SBarry Smith     bases[i] = ldims[i-1];
38147c6ae99SBarry Smith   }
38247c6ae99SBarry Smith   for (i=1; i<=size; i++) {
38347c6ae99SBarry Smith     bases[i] += bases[i-1];
38447c6ae99SBarry Smith   }
385ce00eea3SSatish Balay   base = bases[rank]*dof;
38647c6ae99SBarry Smith 
38747c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
388ce00eea3SSatish Balay   dd->Nlocal = x*y*z*dof;
38947c6ae99SBarry Smith   ierr = VecCreateMPIWithArray(comm,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr);
39047c6ae99SBarry Smith   ierr = VecSetBlockSize(global,dof);CHKERRQ(ierr);
391ce00eea3SSatish Balay   dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs)*dof;
392ce00eea3SSatish Balay   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dd->nlocal,0,&local);CHKERRQ(ierr);
39347c6ae99SBarry Smith   ierr = VecSetBlockSize(local,dof);CHKERRQ(ierr);
39447c6ae99SBarry Smith 
39547c6ae99SBarry Smith   /* generate appropriate vector scatters */
39647c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
39747c6ae99SBarry Smith   ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr);
398ce00eea3SSatish Balay   ierr = ISCreateStride(comm,x*y*z*dof,start,1,&to);CHKERRQ(ierr);
39947c6ae99SBarry Smith 
400db87c5ecSEthan Coon   count = x*y*z;
401ce00eea3SSatish Balay   ierr = PetscMalloc(x*y*z*sizeof(PetscInt),&idx);CHKERRQ(ierr);
402ce00eea3SSatish Balay   left   = xs - Xs; right = left + x;
40347c6ae99SBarry Smith   bottom = ys - Ys; top = bottom + y;
40447c6ae99SBarry Smith   down   = zs - Zs; up  = down + z;
40547c6ae99SBarry Smith   count  = 0;
40647c6ae99SBarry Smith   for (i=down; i<up; i++) {
40747c6ae99SBarry Smith     for (j=bottom; j<top; j++) {
408ce00eea3SSatish Balay       for (k=left; k<right; k++) {
409ce00eea3SSatish Balay         idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
41047c6ae99SBarry Smith       }
41147c6ae99SBarry Smith     }
41247c6ae99SBarry Smith   }
41347c6ae99SBarry Smith 
41447c6ae99SBarry Smith   ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
41547c6ae99SBarry Smith   ierr = VecScatterCreate(local,from,global,to,&ltog);CHKERRQ(ierr);
41647c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,ltog);CHKERRQ(ierr);
417fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
418fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
41947c6ae99SBarry Smith 
420ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
421ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
422aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_BOX) {
423db87c5ecSEthan Coon     count = (IXe-IXs)*(IYe-IYs)*(IZe-IZs);
424db87c5ecSEthan Coon     ierr  = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
425ce00eea3SSatish Balay 
426ce00eea3SSatish Balay     left   = IXs - Xs; right = left + (IXe-IXs);
427ce00eea3SSatish Balay     bottom = IYs - Ys; top = bottom + (IYe-IYs);
428ce00eea3SSatish Balay     down   = IZs - Zs; up  = down + (IZe-IZs);
429ce00eea3SSatish Balay     count = 0;
430ce00eea3SSatish Balay     for (i=down; i<up; i++) {
431ce00eea3SSatish Balay       for (j=bottom; j<top; j++) {
432ce00eea3SSatish Balay         for (k=left; k<right; k++) {
433ce00eea3SSatish Balay           idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
434ce00eea3SSatish Balay         }
435ce00eea3SSatish Balay       }
436ce00eea3SSatish Balay     }
437ce00eea3SSatish Balay     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
438ce00eea3SSatish Balay 
43947c6ae99SBarry Smith   } else {
44047c6ae99SBarry Smith     /* This is way ugly! We need to list the funny cross type region */
441db87c5ecSEthan Coon     count = ((ys-IYs) + (IYe-ye))*x*z + ((xs-IXs) + (IXe-xe))*y*z + ((zs-IZs) + (IZe-ze))*x*y + x*y*z;
442db87c5ecSEthan Coon     ierr   = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
443ce00eea3SSatish Balay 
444ce00eea3SSatish Balay     left   = xs - Xs; right = left + x;
44547c6ae99SBarry Smith     bottom = ys - Ys; top = bottom + y;
44647c6ae99SBarry Smith     down   = zs - Zs;   up  = down + z;
44747c6ae99SBarry Smith     count  = 0;
448ce00eea3SSatish Balay     /* the bottom chunck */
449ce00eea3SSatish Balay     for (i=(IZs-Zs); i<down; i++) {
45047c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
451ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
45247c6ae99SBarry Smith       }
45347c6ae99SBarry Smith     }
45447c6ae99SBarry Smith     /* the middle piece */
45547c6ae99SBarry Smith     for (i=down; i<up; i++) {
45647c6ae99SBarry Smith       /* front */
457ce00eea3SSatish Balay       for (j=(IYs-Ys); j<bottom; j++) {
458ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
45947c6ae99SBarry Smith       }
46047c6ae99SBarry Smith       /* middle */
46147c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
462ce00eea3SSatish Balay         for (k=IXs-Xs; k<IXe-Xs; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
46347c6ae99SBarry Smith       }
46447c6ae99SBarry Smith       /* back */
465ce00eea3SSatish Balay       for (j=top; j<top+IYe-ye; j++) {
466ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
46747c6ae99SBarry Smith       }
46847c6ae99SBarry Smith     }
46947c6ae99SBarry Smith     /* the top piece */
470ce00eea3SSatish Balay     for (i=up; i<up+IZe-ze; 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     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
47647c6ae99SBarry Smith   }
47747c6ae99SBarry Smith 
47847c6ae99SBarry Smith   /* determine who lies on each side of use stored in    n24 n25 n26
47947c6ae99SBarry Smith                                                          n21 n22 n23
48047c6ae99SBarry Smith                                                          n18 n19 n20
48147c6ae99SBarry Smith 
48247c6ae99SBarry Smith                                                          n15 n16 n17
48347c6ae99SBarry Smith                                                          n12     n14
48447c6ae99SBarry Smith                                                          n9  n10 n11
48547c6ae99SBarry Smith 
48647c6ae99SBarry Smith                                                          n6  n7  n8
48747c6ae99SBarry Smith                                                          n3  n4  n5
48847c6ae99SBarry Smith                                                          n0  n1  n2
48947c6ae99SBarry Smith   */
49047c6ae99SBarry Smith 
49147c6ae99SBarry Smith   /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
49247c6ae99SBarry Smith   /* Assume Nodes are Internal to the Cube */
49347c6ae99SBarry Smith   n0  = rank - m*n - m - 1;
49447c6ae99SBarry Smith   n1  = rank - m*n - m;
49547c6ae99SBarry Smith   n2  = rank - m*n - m + 1;
49647c6ae99SBarry Smith   n3  = rank - m*n -1;
49747c6ae99SBarry Smith   n4  = rank - m*n;
49847c6ae99SBarry Smith   n5  = rank - m*n + 1;
49947c6ae99SBarry Smith   n6  = rank - m*n + m - 1;
50047c6ae99SBarry Smith   n7  = rank - m*n + m;
50147c6ae99SBarry Smith   n8  = rank - m*n + m + 1;
50247c6ae99SBarry Smith 
50347c6ae99SBarry Smith   n9  = rank - m - 1;
50447c6ae99SBarry Smith   n10 = rank - m;
50547c6ae99SBarry Smith   n11 = rank - m + 1;
50647c6ae99SBarry Smith   n12 = rank - 1;
50747c6ae99SBarry Smith   n14 = rank + 1;
50847c6ae99SBarry Smith   n15 = rank + m - 1;
50947c6ae99SBarry Smith   n16 = rank + m;
51047c6ae99SBarry Smith   n17 = rank + m + 1;
51147c6ae99SBarry Smith 
51247c6ae99SBarry Smith   n18 = rank + m*n - m - 1;
51347c6ae99SBarry Smith   n19 = rank + m*n - m;
51447c6ae99SBarry Smith   n20 = rank + m*n - m + 1;
51547c6ae99SBarry Smith   n21 = rank + m*n - 1;
51647c6ae99SBarry Smith   n22 = rank + m*n;
51747c6ae99SBarry Smith   n23 = rank + m*n + 1;
51847c6ae99SBarry Smith   n24 = rank + m*n + m - 1;
51947c6ae99SBarry Smith   n25 = rank + m*n + m;
52047c6ae99SBarry Smith   n26 = rank + m*n + m + 1;
52147c6ae99SBarry Smith 
52247c6ae99SBarry Smith   /* Assume Pieces are on Faces of Cube */
52347c6ae99SBarry Smith 
52447c6ae99SBarry Smith   if (xs == 0) { /* First assume not corner or edge */
52547c6ae99SBarry Smith     n0  = rank       -1 - (m*n);
52647c6ae99SBarry Smith     n3  = rank + m   -1 - (m*n);
52747c6ae99SBarry Smith     n6  = rank + 2*m -1 - (m*n);
52847c6ae99SBarry Smith     n9  = rank       -1;
52947c6ae99SBarry Smith     n12 = rank + m   -1;
53047c6ae99SBarry Smith     n15 = rank + 2*m -1;
53147c6ae99SBarry Smith     n18 = rank       -1 + (m*n);
53247c6ae99SBarry Smith     n21 = rank + m   -1 + (m*n);
53347c6ae99SBarry Smith     n24 = rank + 2*m -1 + (m*n);
53447c6ae99SBarry Smith    }
53547c6ae99SBarry Smith 
536ce00eea3SSatish Balay   if (xe == M) { /* First assume not corner or edge */
53747c6ae99SBarry Smith     n2  = rank -2*m +1 - (m*n);
53847c6ae99SBarry Smith     n5  = rank - m  +1 - (m*n);
53947c6ae99SBarry Smith     n8  = rank      +1 - (m*n);
54047c6ae99SBarry Smith     n11 = rank -2*m +1;
54147c6ae99SBarry Smith     n14 = rank - m  +1;
54247c6ae99SBarry Smith     n17 = rank      +1;
54347c6ae99SBarry Smith     n20 = rank -2*m +1 + (m*n);
54447c6ae99SBarry Smith     n23 = rank - m  +1 + (m*n);
54547c6ae99SBarry Smith     n26 = rank      +1 + (m*n);
54647c6ae99SBarry Smith   }
54747c6ae99SBarry Smith 
54847c6ae99SBarry Smith   if (ys==0) { /* First assume not corner or edge */
54947c6ae99SBarry Smith     n0  = rank + m * (n-1) -1 - (m*n);
55047c6ae99SBarry Smith     n1  = rank + m * (n-1)    - (m*n);
55147c6ae99SBarry Smith     n2  = rank + m * (n-1) +1 - (m*n);
55247c6ae99SBarry Smith     n9  = rank + m * (n-1) -1;
55347c6ae99SBarry Smith     n10 = rank + m * (n-1);
55447c6ae99SBarry Smith     n11 = rank + m * (n-1) +1;
55547c6ae99SBarry Smith     n18 = rank + m * (n-1) -1 + (m*n);
55647c6ae99SBarry Smith     n19 = rank + m * (n-1)    + (m*n);
55747c6ae99SBarry Smith     n20 = rank + m * (n-1) +1 + (m*n);
55847c6ae99SBarry Smith   }
55947c6ae99SBarry Smith 
56047c6ae99SBarry Smith   if (ye == N) { /* First assume not corner or edge */
56147c6ae99SBarry Smith     n6  = rank - m * (n-1) -1 - (m*n);
56247c6ae99SBarry Smith     n7  = rank - m * (n-1)    - (m*n);
56347c6ae99SBarry Smith     n8  = rank - m * (n-1) +1 - (m*n);
56447c6ae99SBarry Smith     n15 = rank - m * (n-1) -1;
56547c6ae99SBarry Smith     n16 = rank - m * (n-1);
56647c6ae99SBarry Smith     n17 = rank - m * (n-1) +1;
56747c6ae99SBarry Smith     n24 = rank - m * (n-1) -1 + (m*n);
56847c6ae99SBarry Smith     n25 = rank - m * (n-1)    + (m*n);
56947c6ae99SBarry Smith     n26 = rank - m * (n-1) +1 + (m*n);
57047c6ae99SBarry Smith   }
57147c6ae99SBarry Smith 
57247c6ae99SBarry Smith   if (zs == 0) { /* First assume not corner or edge */
57347c6ae99SBarry Smith     n0 = size - (m*n) + rank - m - 1;
57447c6ae99SBarry Smith     n1 = size - (m*n) + rank - m;
57547c6ae99SBarry Smith     n2 = size - (m*n) + rank - m + 1;
57647c6ae99SBarry Smith     n3 = size - (m*n) + rank - 1;
57747c6ae99SBarry Smith     n4 = size - (m*n) + rank;
57847c6ae99SBarry Smith     n5 = size - (m*n) + rank + 1;
57947c6ae99SBarry Smith     n6 = size - (m*n) + rank + m - 1;
58047c6ae99SBarry Smith     n7 = size - (m*n) + rank + m ;
58147c6ae99SBarry Smith     n8 = size - (m*n) + rank + m + 1;
58247c6ae99SBarry Smith   }
58347c6ae99SBarry Smith 
58447c6ae99SBarry Smith   if (ze == P) { /* First assume not corner or edge */
58547c6ae99SBarry Smith     n18 = (m*n) - (size-rank) - m - 1;
58647c6ae99SBarry Smith     n19 = (m*n) - (size-rank) - m;
58747c6ae99SBarry Smith     n20 = (m*n) - (size-rank) - m + 1;
58847c6ae99SBarry Smith     n21 = (m*n) - (size-rank) - 1;
58947c6ae99SBarry Smith     n22 = (m*n) - (size-rank);
59047c6ae99SBarry Smith     n23 = (m*n) - (size-rank) + 1;
59147c6ae99SBarry Smith     n24 = (m*n) - (size-rank) + m - 1;
59247c6ae99SBarry Smith     n25 = (m*n) - (size-rank) + m;
59347c6ae99SBarry Smith     n26 = (m*n) - (size-rank) + m + 1;
59447c6ae99SBarry Smith   }
59547c6ae99SBarry Smith 
59647c6ae99SBarry Smith   if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */
59747c6ae99SBarry Smith     n0 = size - m*n + rank + m-1 - m;
59847c6ae99SBarry Smith     n3 = size - m*n + rank + m-1;
59947c6ae99SBarry Smith     n6 = size - m*n + rank + m-1 + m;
60047c6ae99SBarry Smith   }
60147c6ae99SBarry Smith 
60247c6ae99SBarry Smith   if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */
60347c6ae99SBarry Smith     n18 = m*n - (size - rank) + m-1 - m;
60447c6ae99SBarry Smith     n21 = m*n - (size - rank) + m-1;
60547c6ae99SBarry Smith     n24 = m*n - (size - rank) + m-1 + m;
60647c6ae99SBarry Smith   }
60747c6ae99SBarry Smith 
60847c6ae99SBarry Smith   if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */
60947c6ae99SBarry Smith     n0  = rank + m*n -1 - m*n;
61047c6ae99SBarry Smith     n9  = rank + m*n -1;
61147c6ae99SBarry Smith     n18 = rank + m*n -1 + m*n;
61247c6ae99SBarry Smith   }
61347c6ae99SBarry Smith 
61447c6ae99SBarry Smith   if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */
61547c6ae99SBarry Smith     n6  = rank - m*(n-1) + m-1 - m*n;
61647c6ae99SBarry Smith     n15 = rank - m*(n-1) + m-1;
61747c6ae99SBarry Smith     n24 = rank - m*(n-1) + m-1 + m*n;
61847c6ae99SBarry Smith   }
61947c6ae99SBarry Smith 
620ce00eea3SSatish Balay   if ((xe==M) && (zs==0)) { /* Assume an edge, not corner */
62147c6ae99SBarry Smith     n2 = size - (m*n-rank) - (m-1) - m;
62247c6ae99SBarry Smith     n5 = size - (m*n-rank) - (m-1);
62347c6ae99SBarry Smith     n8 = size - (m*n-rank) - (m-1) + m;
62447c6ae99SBarry Smith   }
62547c6ae99SBarry Smith 
626ce00eea3SSatish Balay   if ((xe==M) && (ze==P)) { /* Assume an edge, not corner */
62747c6ae99SBarry Smith     n20 = m*n - (size - rank) - (m-1) - m;
62847c6ae99SBarry Smith     n23 = m*n - (size - rank) - (m-1);
62947c6ae99SBarry Smith     n26 = m*n - (size - rank) - (m-1) + m;
63047c6ae99SBarry Smith   }
63147c6ae99SBarry Smith 
632ce00eea3SSatish Balay   if ((xe==M) && (ys==0)) { /* Assume an edge, not corner */
63347c6ae99SBarry Smith     n2  = rank + m*(n-1) - (m-1) - m*n;
63447c6ae99SBarry Smith     n11 = rank + m*(n-1) - (m-1);
63547c6ae99SBarry Smith     n20 = rank + m*(n-1) - (m-1) + m*n;
63647c6ae99SBarry Smith   }
63747c6ae99SBarry Smith 
638ce00eea3SSatish Balay   if ((xe==M) && (ye==N)) { /* Assume an edge, not corner */
63947c6ae99SBarry Smith     n8  = rank - m*n +1 - m*n;
64047c6ae99SBarry Smith     n17 = rank - m*n +1;
64147c6ae99SBarry Smith     n26 = rank - m*n +1 + m*n;
64247c6ae99SBarry Smith   }
64347c6ae99SBarry Smith 
64447c6ae99SBarry Smith   if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */
64547c6ae99SBarry Smith     n0 = size - m + rank -1;
64647c6ae99SBarry Smith     n1 = size - m + rank;
64747c6ae99SBarry Smith     n2 = size - m + rank +1;
64847c6ae99SBarry Smith   }
64947c6ae99SBarry Smith 
65047c6ae99SBarry Smith   if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */
65147c6ae99SBarry Smith     n18 = m*n - (size - rank) + m*(n-1) -1;
65247c6ae99SBarry Smith     n19 = m*n - (size - rank) + m*(n-1);
65347c6ae99SBarry Smith     n20 = m*n - (size - rank) + m*(n-1) +1;
65447c6ae99SBarry Smith   }
65547c6ae99SBarry Smith 
65647c6ae99SBarry Smith   if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */
65747c6ae99SBarry Smith     n6 = size - (m*n-rank) - m * (n-1) -1;
65847c6ae99SBarry Smith     n7 = size - (m*n-rank) - m * (n-1);
65947c6ae99SBarry Smith     n8 = size - (m*n-rank) - m * (n-1) +1;
66047c6ae99SBarry Smith   }
66147c6ae99SBarry Smith 
66247c6ae99SBarry Smith   if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */
66347c6ae99SBarry Smith     n24 = rank - (size-m) -1;
66447c6ae99SBarry Smith     n25 = rank - (size-m);
66547c6ae99SBarry Smith     n26 = rank - (size-m) +1;
66647c6ae99SBarry Smith   }
66747c6ae99SBarry Smith 
66847c6ae99SBarry Smith   /* Check for Corners */
66947c6ae99SBarry Smith   if ((xs==0)   && (ys==0) && (zs==0)) { n0  = size -1;}
67047c6ae99SBarry Smith   if ((xs==0)   && (ys==0) && (ze==P)) { n18 = m*n-1;}
67147c6ae99SBarry Smith   if ((xs==0)   && (ye==N) && (zs==0)) { n6  = (size-1)-m*(n-1);}
67247c6ae99SBarry Smith   if ((xs==0)   && (ye==N) && (ze==P)) { n24 = m-1;}
673ce00eea3SSatish Balay   if ((xe==M) && (ys==0) && (zs==0)) { n2  = size-m;}
674ce00eea3SSatish Balay   if ((xe==M) && (ys==0) && (ze==P)) { n20 = m*n-m;}
675ce00eea3SSatish Balay   if ((xe==M) && (ye==N) && (zs==0)) { n8  = size-m*n;}
676ce00eea3SSatish Balay   if ((xe==M) && (ye==N) && (ze==P)) { n26 = 0;}
67747c6ae99SBarry Smith 
67847c6ae99SBarry Smith   /* Check for when not X,Y, and Z Periodic */
67947c6ae99SBarry Smith 
68047c6ae99SBarry Smith   /* If not X periodic */
6811321219cSEthan Coon   if (bx != DMDA_BOUNDARY_PERIODIC) {
68247c6ae99SBarry Smith     if (xs==0)   {n0  = n3  = n6  = n9  = n12 = n15 = n18 = n21 = n24 = -2;}
683ce00eea3SSatish Balay     if (xe==M) {n2  = n5  = n8  = n11 = n14 = n17 = n20 = n23 = n26 = -2;}
68447c6ae99SBarry Smith   }
68547c6ae99SBarry Smith 
68647c6ae99SBarry Smith   /* If not Y periodic */
6871321219cSEthan Coon   if (by != DMDA_BOUNDARY_PERIODIC) {
68847c6ae99SBarry Smith     if (ys==0)   {n0  = n1  = n2  = n9  = n10 = n11 = n18 = n19 = n20 = -2;}
68947c6ae99SBarry Smith     if (ye==N)   {n6  = n7  = n8  = n15 = n16 = n17 = n24 = n25 = n26 = -2;}
69047c6ae99SBarry Smith   }
69147c6ae99SBarry Smith 
69247c6ae99SBarry Smith   /* If not Z periodic */
6931321219cSEthan Coon   if (bz != DMDA_BOUNDARY_PERIODIC) {
69447c6ae99SBarry Smith     if (zs==0)   {n0  = n1  = n2  = n3  = n4  = n5  = n6  = n7  = n8  = -2;}
69547c6ae99SBarry Smith     if (ze==P)   {n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;}
69647c6ae99SBarry Smith   }
69747c6ae99SBarry Smith 
69847c6ae99SBarry Smith   ierr = PetscMalloc(27*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr);
69947c6ae99SBarry Smith   dd->neighbors[0] = n0;
70047c6ae99SBarry Smith   dd->neighbors[1] = n1;
70147c6ae99SBarry Smith   dd->neighbors[2] = n2;
70247c6ae99SBarry Smith   dd->neighbors[3] = n3;
70347c6ae99SBarry Smith   dd->neighbors[4] = n4;
70447c6ae99SBarry Smith   dd->neighbors[5] = n5;
70547c6ae99SBarry Smith   dd->neighbors[6] = n6;
70647c6ae99SBarry Smith   dd->neighbors[7] = n7;
70747c6ae99SBarry Smith   dd->neighbors[8] = n8;
70847c6ae99SBarry Smith   dd->neighbors[9] = n9;
70947c6ae99SBarry Smith   dd->neighbors[10] = n10;
71047c6ae99SBarry Smith   dd->neighbors[11] = n11;
71147c6ae99SBarry Smith   dd->neighbors[12] = n12;
71247c6ae99SBarry Smith   dd->neighbors[13] = rank;
71347c6ae99SBarry Smith   dd->neighbors[14] = n14;
71447c6ae99SBarry Smith   dd->neighbors[15] = n15;
71547c6ae99SBarry Smith   dd->neighbors[16] = n16;
71647c6ae99SBarry Smith   dd->neighbors[17] = n17;
71747c6ae99SBarry Smith   dd->neighbors[18] = n18;
71847c6ae99SBarry Smith   dd->neighbors[19] = n19;
71947c6ae99SBarry Smith   dd->neighbors[20] = n20;
72047c6ae99SBarry Smith   dd->neighbors[21] = n21;
72147c6ae99SBarry Smith   dd->neighbors[22] = n22;
72247c6ae99SBarry Smith   dd->neighbors[23] = n23;
72347c6ae99SBarry Smith   dd->neighbors[24] = n24;
72447c6ae99SBarry Smith   dd->neighbors[25] = n25;
72547c6ae99SBarry Smith   dd->neighbors[26] = n26;
72647c6ae99SBarry Smith 
72747c6ae99SBarry Smith   /* If star stencil then delete the corner neighbors */
728aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_STAR) {
72947c6ae99SBarry Smith      /* save information about corner neighbors */
73047c6ae99SBarry Smith      sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7;
73147c6ae99SBarry Smith      sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18;
73247c6ae99SBarry Smith      sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25;
73347c6ae99SBarry Smith      sn26 = n26;
73447c6ae99SBarry Smith      n0  = n1  = n2  = n3  = n5  = n6  = n7  = n8  = n9  = n11 =
73547c6ae99SBarry Smith      n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1;
73647c6ae99SBarry Smith   }
73747c6ae99SBarry Smith 
73847c6ae99SBarry Smith 
73947c6ae99SBarry Smith   ierr = PetscMalloc((Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
74047c6ae99SBarry Smith   ierr = PetscLogObjectMemory(da,(Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt));CHKERRQ(ierr);
74147c6ae99SBarry Smith 
74247c6ae99SBarry Smith   nn = 0;
74347c6ae99SBarry Smith   /* Bottom Level */
74447c6ae99SBarry Smith   for (k=0; k<s_z; k++) {
74547c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
74647c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
747ce00eea3SSatish Balay         x_t = lx[n0 % m];
74847c6ae99SBarry Smith         y_t = ly[(n0 % (m*n))/m];
74947c6ae99SBarry Smith         z_t = lz[n0 / (m*n)];
75047c6ae99SBarry 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;
751*6f951b95Secoon         if (twod && (s_t < 0)) {s_t = bases[n0] + x_t*y_t*z_t - (s_y-i)*x_t - s_x;} /* 2D case */
75247c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
75347c6ae99SBarry Smith       }
75447c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
75547c6ae99SBarry Smith         x_t = x;
75647c6ae99SBarry Smith         y_t = ly[(n1 % (m*n))/m];
75747c6ae99SBarry Smith         z_t = lz[n1 / (m*n)];
75847c6ae99SBarry 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;
759*6f951b95Secoon         if (twod && (s_t < 0)) {s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t;} /* 2D case */
76047c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
76147c6ae99SBarry Smith       }
76247c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
763ce00eea3SSatish Balay         x_t = lx[n2 % m];
76447c6ae99SBarry Smith         y_t = ly[(n2 % (m*n))/m];
76547c6ae99SBarry Smith         z_t = lz[n2 / (m*n)];
76647c6ae99SBarry 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;
767*6f951b95Secoon         if (twod && (s_t < 0)) {s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t;} /* 2D case */
76847c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
76947c6ae99SBarry Smith       }
77047c6ae99SBarry Smith     }
77147c6ae99SBarry Smith 
77247c6ae99SBarry Smith     for (i=0; i<y; i++) {
77347c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
774ce00eea3SSatish Balay         x_t = lx[n3 % m];
77547c6ae99SBarry Smith         y_t = y;
77647c6ae99SBarry Smith         z_t = lz[n3 / (m*n)];
77747c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
778*6f951b95Secoon         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 */
77947c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
78047c6ae99SBarry Smith       }
78147c6ae99SBarry Smith 
78247c6ae99SBarry Smith       if (n4 >= 0) { /* middle */
78347c6ae99SBarry Smith         x_t = x;
78447c6ae99SBarry Smith         y_t = y;
78547c6ae99SBarry Smith         z_t = lz[n4 / (m*n)];
78647c6ae99SBarry Smith         s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
787*6f951b95Secoon         if (twod && (s_t < 0)) {s_t = bases[n4] + i*x_t + x_t*y_t*z_t - x_t*y_t;} /* 2D case */
78847c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
78947c6ae99SBarry Smith       }
79047c6ae99SBarry Smith 
79147c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
792ce00eea3SSatish Balay         x_t = lx[n5 % m];
79347c6ae99SBarry Smith         y_t = y;
79447c6ae99SBarry Smith         z_t = lz[n5 / (m*n)];
79547c6ae99SBarry Smith         s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
796*6f951b95Secoon         if (twod && (s_t < 0)) {s_t = bases[n5] + i*x_t + x_t*y_t*z_t - x_t*y_t;} /* 2D case */
79747c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
79847c6ae99SBarry Smith       }
79947c6ae99SBarry Smith     }
80047c6ae99SBarry Smith 
80147c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
80247c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
803ce00eea3SSatish Balay         x_t = lx[n6 % m];
80447c6ae99SBarry Smith         y_t = ly[(n6 % (m*n))/m];
80547c6ae99SBarry Smith         z_t = lz[n6 / (m*n)];
80647c6ae99SBarry Smith         s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
807*6f951b95Secoon         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 */
80847c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
80947c6ae99SBarry Smith       }
81047c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
81147c6ae99SBarry Smith         x_t = x;
81247c6ae99SBarry Smith         y_t = ly[(n7 % (m*n))/m];
81347c6ae99SBarry Smith         z_t = lz[n7 / (m*n)];
81447c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
815*6f951b95Secoon         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 */
81647c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
81747c6ae99SBarry Smith       }
81847c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
819ce00eea3SSatish Balay         x_t = lx[n8 % m];
82047c6ae99SBarry Smith         y_t = ly[(n8 % (m*n))/m];
82147c6ae99SBarry Smith         z_t = lz[n8 / (m*n)];
82247c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
823*6f951b95Secoon         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 */
82447c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
82547c6ae99SBarry Smith       }
82647c6ae99SBarry Smith     }
82747c6ae99SBarry Smith   }
82847c6ae99SBarry Smith 
82947c6ae99SBarry Smith   /* Middle Level */
83047c6ae99SBarry Smith   for (k=0; k<z; k++) {
83147c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
83247c6ae99SBarry Smith       if (n9 >= 0) { /* left below */
833ce00eea3SSatish Balay         x_t = lx[n9 % m];
83447c6ae99SBarry Smith         y_t = ly[(n9 % (m*n))/m];
83547c6ae99SBarry Smith         /* z_t = z; */
83647c6ae99SBarry Smith         s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
83747c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
83847c6ae99SBarry Smith       }
83947c6ae99SBarry Smith       if (n10 >= 0) { /* directly below */
84047c6ae99SBarry Smith         x_t = x;
84147c6ae99SBarry Smith         y_t = ly[(n10 % (m*n))/m];
84247c6ae99SBarry Smith         /* z_t = z; */
84347c6ae99SBarry Smith         s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
84447c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
84547c6ae99SBarry Smith       }
84647c6ae99SBarry Smith       if (n11 >= 0) { /* right below */
847ce00eea3SSatish Balay         x_t = lx[n11 % m];
84847c6ae99SBarry Smith         y_t = ly[(n11 % (m*n))/m];
84947c6ae99SBarry Smith         /* z_t = z; */
85047c6ae99SBarry Smith         s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
85147c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
85247c6ae99SBarry Smith       }
85347c6ae99SBarry Smith     }
85447c6ae99SBarry Smith 
85547c6ae99SBarry Smith     for (i=0; i<y; i++) {
85647c6ae99SBarry Smith       if (n12 >= 0) { /* directly left */
857ce00eea3SSatish Balay         x_t = lx[n12 % m];
85847c6ae99SBarry Smith         y_t = y;
85947c6ae99SBarry Smith         /* z_t = z; */
86047c6ae99SBarry Smith         s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
86147c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
86247c6ae99SBarry Smith       }
86347c6ae99SBarry Smith 
86447c6ae99SBarry Smith       /* Interior */
86547c6ae99SBarry Smith       s_t = bases[rank] + i*x + k*x*y;
86647c6ae99SBarry Smith       for (j=0; j<x; j++) { idx[nn++] = s_t++;}
86747c6ae99SBarry Smith 
86847c6ae99SBarry Smith       if (n14 >= 0) { /* directly right */
869ce00eea3SSatish Balay         x_t = lx[n14 % m];
87047c6ae99SBarry Smith         y_t = y;
87147c6ae99SBarry Smith         /* z_t = z; */
87247c6ae99SBarry Smith         s_t = bases[n14] + i*x_t + k*x_t*y_t;
87347c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
87447c6ae99SBarry Smith       }
87547c6ae99SBarry Smith     }
87647c6ae99SBarry Smith 
87747c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
87847c6ae99SBarry Smith       if (n15 >= 0) { /* left above */
879ce00eea3SSatish Balay         x_t = lx[n15 % m];
88047c6ae99SBarry Smith         y_t = ly[(n15 % (m*n))/m];
88147c6ae99SBarry Smith         /* z_t = z; */
88247c6ae99SBarry Smith         s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
88347c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
88447c6ae99SBarry Smith       }
88547c6ae99SBarry Smith       if (n16 >= 0) { /* directly above */
88647c6ae99SBarry Smith         x_t = x;
88747c6ae99SBarry Smith         y_t = ly[(n16 % (m*n))/m];
88847c6ae99SBarry Smith         /* z_t = z; */
88947c6ae99SBarry Smith         s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
89047c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
89147c6ae99SBarry Smith       }
89247c6ae99SBarry Smith       if (n17 >= 0) { /* right above */
893ce00eea3SSatish Balay         x_t = lx[n17 % m];
89447c6ae99SBarry Smith         y_t = ly[(n17 % (m*n))/m];
89547c6ae99SBarry Smith         /* z_t = z; */
89647c6ae99SBarry Smith         s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
89747c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
89847c6ae99SBarry Smith       }
89947c6ae99SBarry Smith     }
90047c6ae99SBarry Smith   }
90147c6ae99SBarry Smith 
90247c6ae99SBarry Smith   /* Upper Level */
90347c6ae99SBarry Smith   for (k=0; k<s_z; k++) {
90447c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
90547c6ae99SBarry Smith       if (n18 >= 0) { /* left below */
906ce00eea3SSatish Balay         x_t = lx[n18 % m];
90747c6ae99SBarry Smith         y_t = ly[(n18 % (m*n))/m];
90847c6ae99SBarry Smith         /* z_t = lz[n18 / (m*n)]; */
90947c6ae99SBarry Smith         s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
910*6f951b95Secoon         if (twod && (s_t >= M*N*P)) {s_t = bases[n18] - (s_y-i)*x_t -s_x + x_t*y_t;} /* 2d case */
91147c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
91247c6ae99SBarry Smith       }
91347c6ae99SBarry Smith       if (n19 >= 0) { /* directly below */
91447c6ae99SBarry Smith         x_t = x;
91547c6ae99SBarry Smith         y_t = ly[(n19 % (m*n))/m];
91647c6ae99SBarry Smith         /* z_t = lz[n19 / (m*n)]; */
91747c6ae99SBarry Smith         s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
918*6f951b95Secoon         if (twod && (s_t >= M*N*P)) {s_t = bases[n19] - (s_y+1-i)*x_t + x_t*y_t;} /* 2d case */
91947c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
92047c6ae99SBarry Smith       }
92147c6ae99SBarry Smith       if (n20 >= 0) { /* right below */
922ce00eea3SSatish Balay         x_t = lx[n20 % m];
92347c6ae99SBarry Smith         y_t = ly[(n20 % (m*n))/m];
92447c6ae99SBarry Smith         /* z_t = lz[n20 / (m*n)]; */
92547c6ae99SBarry Smith         s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
926*6f951b95Secoon         if (twod && (s_t >= M*N*P)) {s_t = bases[n20] - (s_y+1-i)*x_t + x_t*y_t;} /* 2d case */
92747c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
92847c6ae99SBarry Smith       }
92947c6ae99SBarry Smith     }
93047c6ae99SBarry Smith 
93147c6ae99SBarry Smith     for (i=0; i<y; i++) {
93247c6ae99SBarry Smith       if (n21 >= 0) { /* directly left */
933ce00eea3SSatish Balay         x_t = lx[n21 % m];
93447c6ae99SBarry Smith         y_t = y;
93547c6ae99SBarry Smith         /* z_t = lz[n21 / (m*n)]; */
93647c6ae99SBarry Smith         s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
937*6f951b95Secoon         if (twod && (s_t >= M*N*P)) {s_t = bases[n21] + (i+1)*x_t - s_x;}  /* 2d case */
93847c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
93947c6ae99SBarry Smith       }
94047c6ae99SBarry Smith 
94147c6ae99SBarry Smith       if (n22 >= 0) { /* middle */
94247c6ae99SBarry Smith         x_t = x;
94347c6ae99SBarry Smith         y_t = y;
94447c6ae99SBarry Smith         /* z_t = lz[n22 / (m*n)]; */
94547c6ae99SBarry Smith         s_t = bases[n22] + i*x_t + k*x_t*y_t;
946*6f951b95Secoon         if (twod && (s_t >= M*N*P)) {s_t = bases[n22] + i*x_t;} /* 2d case */
94747c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
94847c6ae99SBarry Smith       }
94947c6ae99SBarry Smith 
95047c6ae99SBarry Smith       if (n23 >= 0) { /* directly right */
951ce00eea3SSatish Balay         x_t = lx[n23 % m];
95247c6ae99SBarry Smith         y_t = y;
95347c6ae99SBarry Smith         /* z_t = lz[n23 / (m*n)]; */
95447c6ae99SBarry Smith         s_t = bases[n23] + i*x_t + k*x_t*y_t;
955*6f951b95Secoon         if (twod && (s_t >= M*N*P)) {s_t = bases[n23] + i*x_t;} /* 2d case */
95647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
95747c6ae99SBarry Smith       }
95847c6ae99SBarry Smith     }
95947c6ae99SBarry Smith 
96047c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
96147c6ae99SBarry Smith       if (n24 >= 0) { /* left above */
962ce00eea3SSatish Balay         x_t = lx[n24 % m];
96347c6ae99SBarry Smith         y_t = ly[(n24 % (m*n))/m];
96447c6ae99SBarry Smith         /* z_t = lz[n24 / (m*n)]; */
96547c6ae99SBarry Smith         s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
966*6f951b95Secoon         if (twod && (s_t >= M*N*P)) {s_t = bases[n24] + i*x_t - s_x;} /* 2d case */
96747c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
96847c6ae99SBarry Smith       }
96947c6ae99SBarry Smith       if (n25 >= 0) { /* directly above */
97047c6ae99SBarry Smith         x_t = x;
97147c6ae99SBarry Smith         y_t = ly[(n25 % (m*n))/m];
97247c6ae99SBarry Smith         /* z_t = lz[n25 / (m*n)]; */
97347c6ae99SBarry Smith         s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
974*6f951b95Secoon         if (twod && (s_t >= M*N*P)) {s_t = bases[n25] + (i-1)*x_t;} /* 2d case */
97547c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
97647c6ae99SBarry Smith       }
97747c6ae99SBarry Smith       if (n26 >= 0) { /* right above */
978ce00eea3SSatish Balay         x_t = lx[n26 % m];
97947c6ae99SBarry Smith         y_t = ly[(n26 % (m*n))/m];
98047c6ae99SBarry Smith         /* z_t = lz[n26 / (m*n)]; */
98147c6ae99SBarry Smith         s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
982*6f951b95Secoon         if (twod && (s_t >= M*N*P)) {s_t = bases[n26] + (i-1)*x_t;} /* 2d case */
98347c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
98447c6ae99SBarry Smith       }
98547c6ae99SBarry Smith     }
98647c6ae99SBarry Smith   }
987ce00eea3SSatish Balay 
988ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
98947c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
99047c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr);
991fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
992fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
99347c6ae99SBarry Smith 
994aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_STAR) {
99547c6ae99SBarry Smith     n0  = sn0;  n1  = sn1;  n2  = sn2;  n3  = sn3;  n5  = sn5;  n6  = sn6; n7 = sn7;
99647c6ae99SBarry Smith     n8  = sn8;  n9  = sn9;  n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18;
99747c6ae99SBarry Smith     n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25;
99847c6ae99SBarry Smith     n26 = sn26;
999ce00eea3SSatish Balay   }
100047c6ae99SBarry Smith 
1001ce00eea3SSatish Balay   if ((stencil_type == DMDA_STENCIL_STAR) ||
10021321219cSEthan Coon       (bx != DMDA_BOUNDARY_PERIODIC && bx) ||
10031321219cSEthan Coon       (by != DMDA_BOUNDARY_PERIODIC && by) ||
10041321219cSEthan Coon       (bz != DMDA_BOUNDARY_PERIODIC && bz)) {
1005ce00eea3SSatish Balay     /*
1006ce00eea3SSatish Balay         Recompute the local to global mappings, this time keeping the
1007ce00eea3SSatish Balay       information about the cross corner processor numbers.
1008ce00eea3SSatish Balay     */
100947c6ae99SBarry Smith     nn = 0;
101047c6ae99SBarry Smith     /* Bottom Level */
101147c6ae99SBarry Smith     for (k=0; k<s_z; k++) {
101247c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
101347c6ae99SBarry Smith         if (n0 >= 0) { /* left below */
1014ce00eea3SSatish Balay           x_t = lx[n0 % m];
101547c6ae99SBarry Smith           y_t = ly[(n0 % (m*n))/m];
101647c6ae99SBarry Smith           z_t = lz[n0 / (m*n)];
101747c6ae99SBarry 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;
101847c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1019ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Ys-ys < 0 && Zs-zs < 0) {
1020ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
102147c6ae99SBarry Smith         }
102247c6ae99SBarry Smith         if (n1 >= 0) { /* directly below */
102347c6ae99SBarry Smith           x_t = x;
102447c6ae99SBarry Smith           y_t = ly[(n1 % (m*n))/m];
102547c6ae99SBarry Smith           z_t = lz[n1 / (m*n)];
102647c6ae99SBarry 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;
102747c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1028ce00eea3SSatish Balay         } else if (Ys-ys < 0 && Zs-zs < 0) {
1029ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
103047c6ae99SBarry Smith         }
103147c6ae99SBarry Smith         if (n2 >= 0) { /* right below */
1032ce00eea3SSatish Balay           x_t = lx[n2 % m];
103347c6ae99SBarry Smith           y_t = ly[(n2 % (m*n))/m];
103447c6ae99SBarry Smith           z_t = lz[n2 / (m*n)];
103547c6ae99SBarry 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;
103647c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1037ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Ys-ys < 0 && Zs-zs < 0) {
1038ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
103947c6ae99SBarry Smith         }
104047c6ae99SBarry Smith       }
104147c6ae99SBarry Smith 
104247c6ae99SBarry Smith       for (i=0; i<y; i++) {
104347c6ae99SBarry Smith         if (n3 >= 0) { /* directly left */
1044ce00eea3SSatish Balay           x_t = lx[n3 % m];
104547c6ae99SBarry Smith           y_t = y;
104647c6ae99SBarry Smith           z_t = lz[n3 / (m*n)];
104747c6ae99SBarry Smith           s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
104847c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1049ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Zs-zs < 0) {
1050ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
105147c6ae99SBarry Smith         }
105247c6ae99SBarry Smith 
105347c6ae99SBarry Smith         if (n4 >= 0) { /* middle */
105447c6ae99SBarry Smith           x_t = x;
105547c6ae99SBarry Smith           y_t = y;
105647c6ae99SBarry Smith           z_t = lz[n4 / (m*n)];
105747c6ae99SBarry Smith           s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
105847c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1059ce00eea3SSatish Balay         } else if (Zs-zs < 0) {
1060ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
106147c6ae99SBarry Smith         }
106247c6ae99SBarry Smith 
106347c6ae99SBarry Smith         if (n5 >= 0) { /* directly right */
1064ce00eea3SSatish Balay           x_t = lx[n5 % m];
106547c6ae99SBarry Smith           y_t = y;
106647c6ae99SBarry Smith           z_t = lz[n5 / (m*n)];
106747c6ae99SBarry Smith           s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
106847c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1069ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Zs-zs < 0) {
1070ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
107147c6ae99SBarry Smith         }
107247c6ae99SBarry Smith       }
107347c6ae99SBarry Smith 
107447c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
107547c6ae99SBarry Smith         if (n6 >= 0) { /* left above */
1076ce00eea3SSatish Balay           x_t = lx[n6 % m];
107747c6ae99SBarry Smith           y_t = ly[(n6 % (m*n))/m];
107847c6ae99SBarry Smith           z_t = lz[n6 / (m*n)];
107947c6ae99SBarry Smith           s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
108047c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1081ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ye-Ye < 0 && Zs-zs < 0) {
1082ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
108347c6ae99SBarry Smith         }
108447c6ae99SBarry Smith         if (n7 >= 0) { /* directly above */
108547c6ae99SBarry Smith           x_t = x;
108647c6ae99SBarry Smith           y_t = ly[(n7 % (m*n))/m];
108747c6ae99SBarry Smith           z_t = lz[n7 / (m*n)];
108847c6ae99SBarry Smith           s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
108947c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1090ce00eea3SSatish Balay         } else if (ye-Ye < 0 && Zs-zs < 0) {
1091ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
109247c6ae99SBarry Smith         }
109347c6ae99SBarry Smith         if (n8 >= 0) { /* right above */
1094ce00eea3SSatish Balay           x_t = lx[n8 % m];
109547c6ae99SBarry Smith           y_t = ly[(n8 % (m*n))/m];
109647c6ae99SBarry Smith           z_t = lz[n8 / (m*n)];
109747c6ae99SBarry Smith           s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
109847c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1099ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ye-Ye < 0 && Zs-zs < 0) {
1100ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
110147c6ae99SBarry Smith         }
110247c6ae99SBarry Smith       }
110347c6ae99SBarry Smith     }
110447c6ae99SBarry Smith 
110547c6ae99SBarry Smith     /* Middle Level */
110647c6ae99SBarry Smith     for (k=0; k<z; k++) {
110747c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
110847c6ae99SBarry Smith         if (n9 >= 0) { /* left below */
1109ce00eea3SSatish Balay           x_t = lx[n9 % m];
111047c6ae99SBarry Smith           y_t = ly[(n9 % (m*n))/m];
111147c6ae99SBarry Smith           /* z_t = z; */
111247c6ae99SBarry Smith           s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
111347c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1114ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Ys-ys < 0) {
1115ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
111647c6ae99SBarry Smith         }
111747c6ae99SBarry Smith         if (n10 >= 0) { /* directly below */
111847c6ae99SBarry Smith           x_t = x;
111947c6ae99SBarry Smith           y_t = ly[(n10 % (m*n))/m];
112047c6ae99SBarry Smith           /* z_t = z; */
112147c6ae99SBarry Smith           s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
112247c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1123ce00eea3SSatish Balay         } else if (Ys-ys < 0) {
1124ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
112547c6ae99SBarry Smith         }
112647c6ae99SBarry Smith         if (n11 >= 0) { /* right below */
1127ce00eea3SSatish Balay           x_t = lx[n11 % m];
112847c6ae99SBarry Smith           y_t = ly[(n11 % (m*n))/m];
112947c6ae99SBarry Smith           /* z_t = z; */
113047c6ae99SBarry Smith           s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
113147c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1132ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Ys-ys < 0) {
1133ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
113447c6ae99SBarry Smith         }
113547c6ae99SBarry Smith       }
113647c6ae99SBarry Smith 
113747c6ae99SBarry Smith       for (i=0; i<y; i++) {
113847c6ae99SBarry Smith         if (n12 >= 0) { /* directly left */
1139ce00eea3SSatish Balay           x_t = lx[n12 % m];
114047c6ae99SBarry Smith           y_t = y;
114147c6ae99SBarry Smith           /* z_t = z; */
114247c6ae99SBarry Smith           s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
114347c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1144ce00eea3SSatish Balay         } else if (Xs-xs < 0) {
1145ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
114647c6ae99SBarry Smith         }
114747c6ae99SBarry Smith 
114847c6ae99SBarry Smith         /* Interior */
114947c6ae99SBarry Smith         s_t = bases[rank] + i*x + k*x*y;
115047c6ae99SBarry Smith         for (j=0; j<x; j++) { idx[nn++] = s_t++;}
115147c6ae99SBarry Smith 
115247c6ae99SBarry Smith         if (n14 >= 0) { /* directly right */
1153ce00eea3SSatish Balay           x_t = lx[n14 % m];
115447c6ae99SBarry Smith           y_t = y;
115547c6ae99SBarry Smith           /* z_t = z; */
115647c6ae99SBarry Smith           s_t = bases[n14] + i*x_t + k*x_t*y_t;
115747c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1158ce00eea3SSatish Balay         } else if (xe-Xe < 0) {
1159ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
116047c6ae99SBarry Smith         }
116147c6ae99SBarry Smith       }
116247c6ae99SBarry Smith 
116347c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
116447c6ae99SBarry Smith         if (n15 >= 0) { /* left above */
1165ce00eea3SSatish Balay           x_t = lx[n15 % m];
116647c6ae99SBarry Smith           y_t = ly[(n15 % (m*n))/m];
116747c6ae99SBarry Smith           /* z_t = z; */
116847c6ae99SBarry Smith           s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
116947c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1170ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ye-Ye < 0) {
1171ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
117247c6ae99SBarry Smith         }
117347c6ae99SBarry Smith         if (n16 >= 0) { /* directly above */
117447c6ae99SBarry Smith           x_t = x;
117547c6ae99SBarry Smith           y_t = ly[(n16 % (m*n))/m];
117647c6ae99SBarry Smith           /* z_t = z; */
117747c6ae99SBarry Smith           s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
117847c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1179ce00eea3SSatish Balay         } else if (ye-Ye < 0) {
1180ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
118147c6ae99SBarry Smith         }
118247c6ae99SBarry Smith         if (n17 >= 0) { /* right above */
1183ce00eea3SSatish Balay           x_t = lx[n17 % m];
118447c6ae99SBarry Smith           y_t = ly[(n17 % (m*n))/m];
118547c6ae99SBarry Smith           /* z_t = z; */
118647c6ae99SBarry Smith           s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
118747c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1188ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ye-Ye < 0) {
1189ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
119047c6ae99SBarry Smith         }
119147c6ae99SBarry Smith       }
119247c6ae99SBarry Smith     }
119347c6ae99SBarry Smith 
119447c6ae99SBarry Smith     /* Upper Level */
119547c6ae99SBarry Smith     for (k=0; k<s_z; k++) {
119647c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
119747c6ae99SBarry Smith         if (n18 >= 0) { /* left below */
1198ce00eea3SSatish Balay           x_t = lx[n18 % m];
119947c6ae99SBarry Smith           y_t = ly[(n18 % (m*n))/m];
120047c6ae99SBarry Smith           /* z_t = lz[n18 / (m*n)]; */
120147c6ae99SBarry Smith           s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
120247c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1203ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Ys-ys < 0 && ze-Ze < 0) {
1204ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
120547c6ae99SBarry Smith         }
120647c6ae99SBarry Smith         if (n19 >= 0) { /* directly below */
120747c6ae99SBarry Smith           x_t = x;
120847c6ae99SBarry Smith           y_t = ly[(n19 % (m*n))/m];
120947c6ae99SBarry Smith           /* z_t = lz[n19 / (m*n)]; */
121047c6ae99SBarry Smith           s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
121147c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1212ce00eea3SSatish Balay         } else if (Ys-ys < 0 && ze-Ze < 0) {
1213ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
121447c6ae99SBarry Smith         }
121547c6ae99SBarry Smith         if (n20 >= 0) { /* right below */
1216ce00eea3SSatish Balay           x_t = lx[n20 % m];
121747c6ae99SBarry Smith           y_t = ly[(n20 % (m*n))/m];
121847c6ae99SBarry Smith           /* z_t = lz[n20 / (m*n)]; */
121947c6ae99SBarry Smith           s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
122047c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1221ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Ys-ys < 0 && ze-Ze < 0) {
1222ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
122347c6ae99SBarry Smith         }
122447c6ae99SBarry Smith       }
122547c6ae99SBarry Smith 
122647c6ae99SBarry Smith       for (i=0; i<y; i++) {
122747c6ae99SBarry Smith         if (n21 >= 0) { /* directly left */
1228ce00eea3SSatish Balay           x_t = lx[n21 % m];
122947c6ae99SBarry Smith           y_t = y;
123047c6ae99SBarry Smith           /* z_t = lz[n21 / (m*n)]; */
123147c6ae99SBarry Smith           s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
123247c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1233ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ze-Ze < 0) {
1234ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
123547c6ae99SBarry Smith         }
123647c6ae99SBarry Smith 
123747c6ae99SBarry Smith         if (n22 >= 0) { /* middle */
123847c6ae99SBarry Smith           x_t = x;
123947c6ae99SBarry Smith           y_t = y;
124047c6ae99SBarry Smith           /* z_t = lz[n22 / (m*n)]; */
124147c6ae99SBarry Smith           s_t = bases[n22] + i*x_t + k*x_t*y_t;
124247c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1243ce00eea3SSatish Balay         } else if (ze-Ze < 0) {
1244ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
124547c6ae99SBarry Smith         }
124647c6ae99SBarry Smith 
124747c6ae99SBarry Smith         if (n23 >= 0) { /* directly right */
1248ce00eea3SSatish Balay           x_t = lx[n23 % m];
124947c6ae99SBarry Smith           y_t = y;
125047c6ae99SBarry Smith           /* z_t = lz[n23 / (m*n)]; */
125147c6ae99SBarry Smith           s_t = bases[n23] + i*x_t + k*x_t*y_t;
125247c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1253ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ze-Ze < 0) {
1254ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
125547c6ae99SBarry Smith         }
125647c6ae99SBarry Smith       }
125747c6ae99SBarry Smith 
125847c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
125947c6ae99SBarry Smith         if (n24 >= 0) { /* left above */
1260ce00eea3SSatish Balay           x_t = lx[n24 % m];
126147c6ae99SBarry Smith           y_t = ly[(n24 % (m*n))/m];
126247c6ae99SBarry Smith           /* z_t = lz[n24 / (m*n)]; */
126347c6ae99SBarry Smith           s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
126447c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1265ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ye-Ye < 0 && ze-Ze < 0) {
1266ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
126747c6ae99SBarry Smith         }
126847c6ae99SBarry Smith         if (n25 >= 0) { /* directly above */
126947c6ae99SBarry Smith           x_t = x;
127047c6ae99SBarry Smith           y_t = ly[(n25 % (m*n))/m];
127147c6ae99SBarry Smith           /* z_t = lz[n25 / (m*n)]; */
127247c6ae99SBarry Smith           s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
127347c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1274ce00eea3SSatish Balay         } else if (ye-Ye < 0 && ze-Ze < 0) {
1275ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
127647c6ae99SBarry Smith         }
127747c6ae99SBarry Smith         if (n26 >= 0) { /* right above */
1278ce00eea3SSatish Balay           x_t = lx[n26 % m];
127947c6ae99SBarry Smith           y_t = ly[(n26 % (m*n))/m];
128047c6ae99SBarry Smith           /* z_t = lz[n26 / (m*n)]; */
128147c6ae99SBarry Smith           s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
128247c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1283ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ye-Ye < 0 && ze-Ze < 0) {
1284ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
128547c6ae99SBarry Smith         }
128647c6ae99SBarry Smith       }
128747c6ae99SBarry Smith     }
128847c6ae99SBarry Smith   }
128947c6ae99SBarry Smith   /*
129047c6ae99SBarry Smith      Set the local to global ordering in the global vector, this allows use
129147c6ae99SBarry Smith      of VecSetValuesLocal().
129247c6ae99SBarry Smith   */
1293ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_OWN_POINTER,&ltogis);CHKERRQ(ierr);
1294ce00eea3SSatish Balay   ierr = PetscMalloc(nn*dof*sizeof(PetscInt),&idx_cpy);CHKERRQ(ierr);
1295db87c5ecSEthan Coon   ierr = PetscLogObjectMemory(da,nn*dof*sizeof(PetscInt));CHKERRQ(ierr);
1296ce00eea3SSatish Balay   ierr = ISGetIndices(ltogis, &idx_full);
1297ce00eea3SSatish Balay   ierr = PetscMemcpy(idx_cpy,idx_full,nn*dof*sizeof(PetscInt));CHKERRQ(ierr);
1298ce00eea3SSatish Balay   ierr = ISRestoreIndices(ltogis, &idx_full);
1299ce00eea3SSatish Balay   ierr = ISLocalToGlobalMappingCreateIS(ltogis,&da->ltogmap);CHKERRQ(ierr);
1300ce00eea3SSatish Balay   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
1301fcfd50ebSBarry Smith   ierr = ISDestroy(&ltogis);CHKERRQ(ierr);
13021411c6eeSJed Brown   ierr = ISLocalToGlobalMappingBlock(da->ltogmap,dd->w,&da->ltogmapb);CHKERRQ(ierr);
13031411c6eeSJed Brown   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
130447c6ae99SBarry Smith 
1305ce00eea3SSatish Balay   ierr = PetscFree2(bases,ldims);CHKERRQ(ierr);
1306ce00eea3SSatish Balay   dd->m  = m;  dd->n  = n;  dd->p  = p;
1307ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
1308ce00eea3SSatish Balay   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = zs; dd->ze = ze;
1309ce00eea3SSatish Balay   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = Zs; dd->Ze = Ze;
1310ce00eea3SSatish Balay 
1311fcfd50ebSBarry Smith   ierr = VecDestroy(&local);CHKERRQ(ierr);
1312fcfd50ebSBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
1313ce00eea3SSatish Balay 
1314ce00eea3SSatish Balay   dd->gtol      = gtol;
1315ce00eea3SSatish Balay   dd->ltog      = ltog;
1316ce00eea3SSatish Balay   dd->idx       = idx_cpy;
1317ce00eea3SSatish Balay   dd->Nl        = nn*dof;
1318ce00eea3SSatish Balay   dd->base      = base;
1319ce00eea3SSatish Balay   da->ops->view = DMView_DA_3d;
132047c6ae99SBarry Smith   dd->ltol = PETSC_NULL;
132147c6ae99SBarry Smith   dd->ao   = PETSC_NULL;
1322ce00eea3SSatish Balay 
132347c6ae99SBarry Smith   PetscFunctionReturn(0);
132447c6ae99SBarry Smith }
1325564755cdSBarry Smith 
132647c6ae99SBarry Smith 
132747c6ae99SBarry Smith #undef __FUNCT__
1328aa219208SBarry Smith #define __FUNCT__ "DMDACreate3d"
132947c6ae99SBarry Smith /*@C
1330aa219208SBarry Smith    DMDACreate3d - Creates an object that will manage the communication of three-dimensional
133147c6ae99SBarry Smith    regular array data that is distributed across some processors.
133247c6ae99SBarry Smith 
133347c6ae99SBarry Smith    Collective on MPI_Comm
133447c6ae99SBarry Smith 
133547c6ae99SBarry Smith    Input Parameters:
133647c6ae99SBarry Smith +  comm - MPI communicator
13371321219cSEthan Coon .  bx,by,bz - type of ghost nodes the array have.
13381321219cSEthan Coon          Use one of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC.
1339aa219208SBarry Smith .  stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX)
134047c6ae99SBarry 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
134147c6ae99SBarry Smith             from the command line with -da_grid_x <M> -da_grid_y <N> -da_grid_z <P>)
134247c6ae99SBarry Smith .  m,n,p - corresponding number of processors in each dimension
134347c6ae99SBarry Smith            (or PETSC_DECIDE to have calculated)
134447c6ae99SBarry Smith .  dof - number of degrees of freedom per node
134547c6ae99SBarry Smith .  lx, ly, lz - arrays containing the number of nodes in each cell along
134647c6ae99SBarry Smith           the x, y, and z coordinates, or PETSC_NULL. If non-null, these
134747c6ae99SBarry Smith           must be of length as m,n,p and the corresponding
134847c6ae99SBarry Smith           m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of
134947c6ae99SBarry Smith           the ly[] must N, sum of the lz[] must be P
135047c6ae99SBarry Smith -  s - stencil width
135147c6ae99SBarry Smith 
135247c6ae99SBarry Smith    Output Parameter:
135347c6ae99SBarry Smith .  da - the resulting distributed array object
135447c6ae99SBarry Smith 
135547c6ae99SBarry Smith    Options Database Key:
1356aa219208SBarry Smith +  -da_view - Calls DMView() at the conclusion of DMDACreate3d()
135747c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
135847c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
135947c6ae99SBarry Smith .  -da_grid_z <nz> - number of grid points in z direction, if P < 0
136047c6ae99SBarry Smith .  -da_processors_x <MX> - number of processors in x direction
136147c6ae99SBarry Smith .  -da_processors_y <MY> - number of processors in y direction
136247c6ae99SBarry Smith .  -da_processors_z <MZ> - number of processors in z direction
1363e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement ratio in x direction
1364e0f5d30fSBarry Smith .  -da_refine_y <ry> - refinement ratio in y direction
1365e0f5d30fSBarry Smith .  -da_refine_z <rz>- refinement ratio in z directio
1366e0f5d30fSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating it, , if M, N, or P < 0
136747c6ae99SBarry Smith 
136847c6ae99SBarry Smith    Level: beginner
136947c6ae99SBarry Smith 
137047c6ae99SBarry Smith    Notes:
1371aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
1372aa219208SBarry Smith    standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
137347c6ae99SBarry Smith    the standard 27-pt stencil.
137447c6ae99SBarry Smith 
1375aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
1376564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
1377564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
137847c6ae99SBarry Smith 
137947c6ae99SBarry Smith .keywords: distributed array, create, three-dimensional
138047c6ae99SBarry Smith 
1381aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate2d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
1382aa219208SBarry Smith           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDASetRefinementFactor(),
1383d461ba97SBarry Smith           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
138447c6ae99SBarry Smith 
138547c6ae99SBarry Smith @*/
13861321219cSEthan Coon PetscErrorCode  DMDACreate3d(MPI_Comm comm,DMDABoundaryType bx,DMDABoundaryType by,DMDABoundaryType bz,DMDAStencilType stencil_type,PetscInt M,
13879a42bb27SBarry 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)
138847c6ae99SBarry Smith {
138947c6ae99SBarry Smith   PetscErrorCode ierr;
139047c6ae99SBarry Smith 
139147c6ae99SBarry Smith   PetscFunctionBegin;
1392aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
1393aa219208SBarry Smith   ierr = DMDASetDim(*da, 3);CHKERRQ(ierr);
1394aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, P);CHKERRQ(ierr);
1395aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, p);CHKERRQ(ierr);
13961321219cSEthan Coon   ierr = DMDASetBoundaryType(*da, bx, by, bz);CHKERRQ(ierr);
1397aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
1398aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
1399aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
1400aa219208SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, lz);CHKERRQ(ierr);
140147c6ae99SBarry Smith   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
14029a42bb27SBarry Smith   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
14039a42bb27SBarry Smith   ierr = DMSetUp(*da);CHKERRQ(ierr);
14047242296bSJed Brown   ierr = DMView_DA_Private(*da);CHKERRQ(ierr);
140547c6ae99SBarry Smith   PetscFunctionReturn(0);
140647c6ae99SBarry Smith }
1407