xref: /petsc/src/dm/impls/da/da3.c (revision 3bf036e263fd78807e2931ff42d9ddcd8aae3fd4)
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 
7b45d2f2cSJed Brown #include <petsc-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 
24251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
25251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
26251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
279a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
28251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((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)
426636e97aSMatthew G Knepley       if (da->coordinates) {
4347c6ae99SBarry Smith         PetscInt        last;
4447c6ae99SBarry Smith         const PetscReal *coors;
456636e97aSMatthew G Knepley         ierr = VecGetArrayRead(da->coordinates,&coors);CHKERRQ(ierr);
466636e97aSMatthew G Knepley         ierr = VecGetLocalSize(da->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);
496636e97aSMatthew G Knepley         ierr = VecRestoreArrayRead(da->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;
17688661749SPeter Brune   PetscInt         o            = dd->overlap;
17747c6ae99SBarry Smith   const PetscInt   dof          = dd->w;
17847c6ae99SBarry Smith   const PetscInt   s            = dd->s;
17919fd82e9SBarry Smith   DMDABoundaryType bx         = dd->bx;
18019fd82e9SBarry Smith   DMDABoundaryType by         = dd->by;
18119fd82e9SBarry Smith   DMDABoundaryType bz         = dd->bz;
18219fd82e9SBarry Smith   DMDAStencilType  stencil_type = dd->stencil_type;
18347c6ae99SBarry Smith   PetscInt         *lx           = dd->lx;
18447c6ae99SBarry Smith   PetscInt         *ly           = dd->ly;
18547c6ae99SBarry Smith   PetscInt         *lz           = dd->lz;
18647c6ae99SBarry Smith   MPI_Comm         comm;
18747c6ae99SBarry Smith   PetscMPIInt      rank,size;
188ce00eea3SSatish Balay   PetscInt         xs = 0,xe,ys = 0,ye,zs = 0,ze,x = 0,y = 0,z = 0;
189ce00eea3SSatish Balay   PetscInt         Xs,Xe,Ys,Ye,Zs,Ze,IXs,IXe,IYs,IYe,IZs,IZe,start,end,pm;
190ce00eea3SSatish Balay   PetscInt         left,right,up,down,bottom,top,i,j,k,*idx,*idx_cpy,nn;
191ce00eea3SSatish Balay   const PetscInt   *idx_full;
19247c6ae99SBarry Smith   PetscInt         n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14;
19347c6ae99SBarry Smith   PetscInt         n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26;
194db87c5ecSEthan Coon   PetscInt         *bases,*ldims,base,x_t,y_t,z_t,s_t,count,s_x,s_y,s_z;
19547c6ae99SBarry Smith   PetscInt         sn0 = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0;
19647c6ae99SBarry Smith   PetscInt         sn8 = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0;
19747c6ae99SBarry Smith   PetscInt         sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0;
19847c6ae99SBarry Smith   Vec              local,global;
19947c6ae99SBarry Smith   VecScatter       ltog,gtol;
200ce00eea3SSatish Balay   IS               to,from,ltogis;
2016f951b95Secoon   PetscBool        twod;
20247c6ae99SBarry Smith   PetscErrorCode   ierr;
20347c6ae99SBarry Smith 
2046f951b95Secoon 
20547c6ae99SBarry Smith   PetscFunctionBegin;
206c2859e5eSBarry Smith   if (stencil_type == DMDA_STENCIL_BOX && (bx == DMDA_BOUNDARY_MIRROR || by == DMDA_BOUNDARY_MIRROR || bz == DMDA_BOUNDARY_MIRROR)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Mirror boundary and box stencil");
20747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr);
2083855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES)
2093855c12bSBarry Smith   if (((Petsc64bitInt) M)*((Petsc64bitInt) N)*((Petsc64bitInt) P)*((Petsc64bitInt) dof) > (Petsc64bitInt) 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);
2103855c12bSBarry Smith #endif
2113855c12bSBarry Smith 
21247c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
21347c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
21447c6ae99SBarry Smith 
21547c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
21647c6ae99SBarry Smith     if (m < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
21747c6ae99SBarry Smith     else if (m > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
21847c6ae99SBarry Smith   }
21947c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
22047c6ae99SBarry Smith     if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
22147c6ae99SBarry Smith     else if (n > size)  SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
22247c6ae99SBarry Smith   }
22347c6ae99SBarry Smith   if (p != PETSC_DECIDE) {
22447c6ae99SBarry Smith     if (p < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %D",p);
22547c6ae99SBarry Smith     else if (p > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %D %d",p,size);
22647c6ae99SBarry Smith   }
2270716a85fSBarry 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);
22847c6ae99SBarry Smith 
22947c6ae99SBarry Smith   /* Partition the array among the processors */
23047c6ae99SBarry Smith   if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) {
23147c6ae99SBarry Smith     m = size/(n*p);
23247c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
23347c6ae99SBarry Smith     n = size/(m*p);
23447c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
23547c6ae99SBarry Smith     p = size/(m*n);
23647c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
23747c6ae99SBarry Smith     /* try for squarish distribution */
2388f1a2a5eSBarry Smith     m = (int)(0.5 + sqrt(((double)M)*((double)size)/((double)N*p)));
23947c6ae99SBarry Smith     if (!m) m = 1;
24047c6ae99SBarry Smith     while (m > 0) {
24147c6ae99SBarry Smith       n = size/(m*p);
24247c6ae99SBarry Smith       if (m*n*p == size) break;
24347c6ae99SBarry Smith       m--;
24447c6ae99SBarry Smith     }
24547c6ae99SBarry Smith     if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %D",p);
24647c6ae99SBarry Smith     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
24747c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
24847c6ae99SBarry Smith     /* try for squarish distribution */
2498f1a2a5eSBarry Smith     m = (int)(0.5 + sqrt(((double)M)*((double)size)/((double)P*n)));
25047c6ae99SBarry Smith     if (!m) m = 1;
25147c6ae99SBarry Smith     while (m > 0) {
25247c6ae99SBarry Smith       p = size/(m*n);
25347c6ae99SBarry Smith       if (m*n*p == size) break;
25447c6ae99SBarry Smith       m--;
25547c6ae99SBarry Smith     }
25647c6ae99SBarry Smith     if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %D",n);
25747c6ae99SBarry Smith     if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;}
25847c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
25947c6ae99SBarry Smith     /* try for squarish distribution */
2608f1a2a5eSBarry Smith     n = (int)(0.5 + sqrt(((double)N)*((double)size)/((double)P*m)));
26147c6ae99SBarry Smith     if (!n) n = 1;
26247c6ae99SBarry Smith     while (n > 0) {
26347c6ae99SBarry Smith       p = size/(m*n);
26447c6ae99SBarry Smith       if (m*n*p == size) break;
26547c6ae99SBarry Smith       n--;
26647c6ae99SBarry Smith     }
26747c6ae99SBarry Smith     if (!n) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %D",n);
26847c6ae99SBarry Smith     if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;}
26947c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
27047c6ae99SBarry Smith     /* try for squarish distribution */
2718f1a2a5eSBarry Smith     n = (PetscInt)(0.5 + pow(((double)N*N)*((double)size)/((double)P*M),(double)(1./3.)));
27247c6ae99SBarry Smith     if (!n) n = 1;
27347c6ae99SBarry Smith     while (n > 0) {
27447c6ae99SBarry Smith       pm = size/n;
27547c6ae99SBarry Smith       if (n*pm == size) break;
27647c6ae99SBarry Smith       n--;
27747c6ae99SBarry Smith     }
27847c6ae99SBarry Smith     if (!n) n = 1;
2798f1a2a5eSBarry Smith     m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)P*n)));
28047c6ae99SBarry Smith     if (!m) m = 1;
28147c6ae99SBarry Smith     while (m > 0) {
28247c6ae99SBarry Smith       p = size/(m*n);
28347c6ae99SBarry Smith       if (m*n*p == size) break;
28447c6ae99SBarry Smith       m--;
28547c6ae99SBarry Smith     }
28647c6ae99SBarry Smith     if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;}
28747c6ae99SBarry Smith   } else if (m*n*p != size) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
28847c6ae99SBarry Smith 
28947c6ae99SBarry Smith   if (m*n*p != size) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_PLIB,"Could not find good partition");
29047c6ae99SBarry Smith   if (M < m) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
29147c6ae99SBarry Smith   if (N < n) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
29247c6ae99SBarry Smith   if (P < p) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %D %D",P,p);
29347c6ae99SBarry Smith 
29447c6ae99SBarry Smith   /*
29547c6ae99SBarry Smith      Determine locally owned region
29647c6ae99SBarry Smith      [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes
29747c6ae99SBarry Smith   */
29847c6ae99SBarry Smith 
29947c6ae99SBarry Smith   if (!lx) {
30047c6ae99SBarry Smith     ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
30147c6ae99SBarry Smith     lx = dd->lx;
30247c6ae99SBarry Smith     for (i=0; i<m; i++) {
30347c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > (i % m));
30447c6ae99SBarry Smith     }
30547c6ae99SBarry Smith   }
30647c6ae99SBarry Smith   x  = lx[rank % m];
30747c6ae99SBarry Smith   xs = 0;
30847c6ae99SBarry Smith   for (i=0; i<(rank%m); i++) { xs += lx[i];}
30988661749SPeter Brune   if ((x < s+o) && ((m > 1) || (bx == DMDA_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);
31047c6ae99SBarry Smith 
31147c6ae99SBarry Smith   if (!ly) {
31247c6ae99SBarry Smith     ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
31347c6ae99SBarry Smith     ly = dd->ly;
31447c6ae99SBarry Smith     for (i=0; i<n; i++) {
31547c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > (i % n));
31647c6ae99SBarry Smith     }
31747c6ae99SBarry Smith   }
31847c6ae99SBarry Smith   y  = ly[(rank % (m*n))/m];
31988661749SPeter Brune   if ((y < s+o) && ((n > 1) || (by == DMDA_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);
32030729d88SBarry Smith 
32147c6ae99SBarry Smith   ys = 0;
32247c6ae99SBarry Smith   for (i=0; i<(rank % (m*n))/m; i++) { ys += ly[i];}
32347c6ae99SBarry Smith 
32447c6ae99SBarry Smith   if (!lz) {
32547c6ae99SBarry Smith     ierr = PetscMalloc(p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr);
32647c6ae99SBarry Smith     lz = dd->lz;
32747c6ae99SBarry Smith     for (i=0; i<p; i++) {
32847c6ae99SBarry Smith       lz[i] = P/p + ((P % p) > (i % p));
32947c6ae99SBarry Smith     }
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
334fdc81ce8SEthan Coon    case when p=P=1 and DMDA_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;
3376f951b95Secoon   if (P == 1) {
3386f951b95Secoon     twod = PETSC_TRUE;
33988661749SPeter Brune   } else if ((z < s+o) && ((p > 1) || (bz == DMDA_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+o);
34047c6ae99SBarry Smith   zs = 0;
34147c6ae99SBarry Smith   for (i=0; i<(rank/(m*n)); i++) { zs += lz[i];}
34247c6ae99SBarry Smith   ye = ys + y;
34347c6ae99SBarry Smith   xe = xs + x;
34447c6ae99SBarry Smith   ze = zs + z;
34547c6ae99SBarry Smith 
34688661749SPeter Brune     /* determine ghost region (Xs) and region scattered into (IXs)  */
34788661749SPeter Brune   if (xs-s-o > 0) {
34888661749SPeter Brune     Xs = xs - s - o; IXs = xs - s - o;
34988661749SPeter Brune   } else {
35088661749SPeter Brune     if (bx) {
35188661749SPeter Brune       Xs = xs - s;
35288661749SPeter Brune     } else {
35388661749SPeter Brune       Xs = 0;
35488661749SPeter Brune     }
35588661749SPeter Brune     IXs = 0;
35688661749SPeter Brune   }
35788661749SPeter Brune   if (xe+s+o <= M) {
35888661749SPeter Brune     Xe = xe + s + o; IXe = xe + s + o;
35988661749SPeter Brune   } else {
36088661749SPeter Brune     if (bx) {
36188661749SPeter Brune       Xs = xs - s - o; Xe = xe + s;
36288661749SPeter Brune     } else {
36388661749SPeter Brune       Xe = M;
36488661749SPeter Brune     }
36588661749SPeter Brune     IXe = M;
36688661749SPeter Brune   }
36747c6ae99SBarry Smith 
368c2859e5eSBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC || bx == DMDA_BOUNDARY_MIRROR) {
36988661749SPeter Brune     IXs = xs - s - o;
37088661749SPeter Brune     IXe = xe + s + o;
37188661749SPeter Brune     Xs = xs - s - o;
37288661749SPeter Brune     Xe = xe + s + o;
37388661749SPeter Brune   }
37488661749SPeter Brune 
37588661749SPeter Brune   if (ys-s-o > 0) {
37688661749SPeter Brune     Ys = ys - s - o; IYs = ys - s - o;
37788661749SPeter Brune   } else {
37888661749SPeter Brune     if (by) {
37988661749SPeter Brune       Ys = ys - s;
38088661749SPeter Brune     } else {
38188661749SPeter Brune       Ys = 0;
38288661749SPeter Brune     }
38388661749SPeter Brune     IYs = 0;
38488661749SPeter Brune   }
38588661749SPeter Brune   if (ye+s+o <= N) {
38688661749SPeter Brune     Ye = ye + s + o; IYe = ye + s + o;
38788661749SPeter Brune   } else {
38888661749SPeter Brune     if (by) {
38988661749SPeter Brune       Ye = ye + s;
39088661749SPeter Brune     } else {
39188661749SPeter Brune       Ye = N;
39288661749SPeter Brune     }
39388661749SPeter Brune     IYe = N;
39488661749SPeter Brune   }
39588661749SPeter Brune 
396c2859e5eSBarry Smith   if (by == DMDA_BOUNDARY_PERIODIC || by == DMDA_BOUNDARY_MIRROR) {
39788661749SPeter Brune     IYs = ys - s - o;
39888661749SPeter Brune     IYe = ye + s + o;
39988661749SPeter Brune     Ys = ys - s - o;
40088661749SPeter Brune     Ye = ye + s + o;
40188661749SPeter Brune   }
40288661749SPeter Brune 
40388661749SPeter Brune   if (zs-s-o > 0) {
40488661749SPeter Brune     Zs = zs - s - o; IZs = zs - s - o;
40588661749SPeter Brune   } else {
40688661749SPeter Brune     if (bz) {
40788661749SPeter Brune       Zs = zs - s;
40888661749SPeter Brune     } else {
40988661749SPeter Brune       Zs = 0;
41088661749SPeter Brune     }
41188661749SPeter Brune     IZs = 0;
41288661749SPeter Brune   }
41388661749SPeter Brune   if (ze+s+o <= P) {
41488661749SPeter Brune     Ze = ze + s + o; IZe = ze + s + o;
41588661749SPeter Brune   } else {
41688661749SPeter Brune     if (bz) {
41788661749SPeter Brune       Ze = ze + s;
41888661749SPeter Brune     } else {
41988661749SPeter Brune       Ze = P;
42088661749SPeter Brune     }
42188661749SPeter Brune     IZe = P;
42288661749SPeter Brune   }
42388661749SPeter Brune 
424c2859e5eSBarry Smith   if (bz == DMDA_BOUNDARY_PERIODIC || bz == DMDA_BOUNDARY_MIRROR) {
42588661749SPeter Brune     IZs = zs - s - o;
42688661749SPeter Brune     IZe = ze + s + o;
42788661749SPeter Brune     Zs = zs - s - o;
42888661749SPeter Brune     Ze = ze + s + o;
42988661749SPeter Brune   }
43047c6ae99SBarry Smith 
43147c6ae99SBarry Smith   /* Resize all X parameters to reflect w */
43288661749SPeter Brune   s_x  = s + o;
43388661749SPeter Brune   s_y  = s + o;
43488661749SPeter Brune   s_z  = s + o;
43547c6ae99SBarry Smith 
43647c6ae99SBarry Smith   /* determine starting point of each processor */
43747c6ae99SBarry Smith   nn       = x*y*z;
43847c6ae99SBarry Smith   ierr     = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr);
43947c6ae99SBarry Smith   ierr     = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
44047c6ae99SBarry Smith   bases[0] = 0;
44147c6ae99SBarry Smith   for (i=1; i<=size; i++) {
44247c6ae99SBarry Smith     bases[i] = ldims[i-1];
44347c6ae99SBarry Smith   }
44447c6ae99SBarry Smith   for (i=1; i<=size; i++) {
44547c6ae99SBarry Smith     bases[i] += bases[i-1];
44647c6ae99SBarry Smith   }
447ce00eea3SSatish Balay   base = bases[rank]*dof;
44847c6ae99SBarry Smith 
44947c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
450ce00eea3SSatish Balay   dd->Nlocal = x*y*z*dof;
451778a2246SBarry Smith   ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr);
452ce00eea3SSatish Balay   dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs)*dof;
453778a2246SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,0,&local);CHKERRQ(ierr);
45447c6ae99SBarry Smith 
45547c6ae99SBarry Smith   /* generate appropriate vector scatters */
45647c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
45747c6ae99SBarry Smith   ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr);
458ce00eea3SSatish Balay   ierr = ISCreateStride(comm,x*y*z*dof,start,1,&to);CHKERRQ(ierr);
45947c6ae99SBarry Smith 
460ce00eea3SSatish Balay   ierr = PetscMalloc(x*y*z*sizeof(PetscInt),&idx);CHKERRQ(ierr);
461ce00eea3SSatish Balay   left   = xs - Xs; right = left + x;
46247c6ae99SBarry Smith   bottom = ys - Ys; top = bottom + y;
46347c6ae99SBarry Smith   down   = zs - Zs; up  = down + z;
46447c6ae99SBarry Smith   count  = 0;
46547c6ae99SBarry Smith   for (i=down; i<up; i++) {
46647c6ae99SBarry Smith     for (j=bottom; j<top; j++) {
467ce00eea3SSatish Balay       for (k=left; k<right; k++) {
468ce00eea3SSatish Balay         idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
46947c6ae99SBarry Smith       }
47047c6ae99SBarry Smith     }
47147c6ae99SBarry Smith   }
47247c6ae99SBarry Smith 
47347c6ae99SBarry Smith   ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
47447c6ae99SBarry Smith   ierr = VecScatterCreate(local,from,global,to,&ltog);CHKERRQ(ierr);
47547c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,ltog);CHKERRQ(ierr);
476fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
477fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
47847c6ae99SBarry Smith 
479ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
480ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
48188661749SPeter Brune   if (stencil_type == DMDA_STENCIL_BOX || o > 0) {
482db87c5ecSEthan Coon     count = (IXe-IXs)*(IYe-IYs)*(IZe-IZs);
483db87c5ecSEthan Coon     ierr  = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
484ce00eea3SSatish Balay 
485ce00eea3SSatish Balay     left   = IXs - Xs; right = left + (IXe-IXs);
486ce00eea3SSatish Balay     bottom = IYs - Ys; top = bottom + (IYe-IYs);
487ce00eea3SSatish Balay     down   = IZs - Zs; up  = down + (IZe-IZs);
488ce00eea3SSatish Balay     count = 0;
489ce00eea3SSatish Balay     for (i=down; i<up; i++) {
490ce00eea3SSatish Balay       for (j=bottom; j<top; j++) {
491ce00eea3SSatish Balay         for (k=left; k<right; k++) {
492ce00eea3SSatish Balay           idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
493ce00eea3SSatish Balay         }
494ce00eea3SSatish Balay       }
495ce00eea3SSatish Balay     }
496ce00eea3SSatish Balay     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
497ce00eea3SSatish Balay 
49847c6ae99SBarry Smith   } else {
49947c6ae99SBarry Smith     /* This is way ugly! We need to list the funny cross type region */
500db87c5ecSEthan Coon     count = ((ys-IYs) + (IYe-ye))*x*z + ((xs-IXs) + (IXe-xe))*y*z + ((zs-IZs) + (IZe-ze))*x*y + x*y*z;
501db87c5ecSEthan Coon     ierr   = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
502ce00eea3SSatish Balay 
503ce00eea3SSatish Balay     left   = xs - Xs; right = left + x;
50447c6ae99SBarry Smith     bottom = ys - Ys; top = bottom + y;
50547c6ae99SBarry Smith     down   = zs - Zs;   up  = down + z;
50647c6ae99SBarry Smith     count  = 0;
507ce00eea3SSatish Balay     /* the bottom chunck */
508ce00eea3SSatish Balay     for (i=(IZs-Zs); i<down; i++) {
50947c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
510ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
51147c6ae99SBarry Smith       }
51247c6ae99SBarry Smith     }
51347c6ae99SBarry Smith     /* the middle piece */
51447c6ae99SBarry Smith     for (i=down; i<up; i++) {
51547c6ae99SBarry Smith       /* front */
516ce00eea3SSatish Balay       for (j=(IYs-Ys); j<bottom; j++) {
517ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
51847c6ae99SBarry Smith       }
51947c6ae99SBarry Smith       /* middle */
52047c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
521ce00eea3SSatish Balay         for (k=IXs-Xs; k<IXe-Xs; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
52247c6ae99SBarry Smith       }
52347c6ae99SBarry Smith       /* back */
524ce00eea3SSatish Balay       for (j=top; j<top+IYe-ye; j++) {
525ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
52647c6ae99SBarry Smith       }
52747c6ae99SBarry Smith     }
52847c6ae99SBarry Smith     /* the top piece */
529ce00eea3SSatish Balay     for (i=up; i<up+IZe-ze; i++) {
53047c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
531ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
53247c6ae99SBarry Smith       }
53347c6ae99SBarry Smith     }
53447c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
53547c6ae99SBarry Smith   }
53647c6ae99SBarry Smith 
53747c6ae99SBarry Smith   /* determine who lies on each side of use stored in    n24 n25 n26
53847c6ae99SBarry Smith                                                          n21 n22 n23
53947c6ae99SBarry Smith                                                          n18 n19 n20
54047c6ae99SBarry Smith 
54147c6ae99SBarry Smith                                                          n15 n16 n17
54247c6ae99SBarry Smith                                                          n12     n14
54347c6ae99SBarry Smith                                                          n9  n10 n11
54447c6ae99SBarry Smith 
54547c6ae99SBarry Smith                                                          n6  n7  n8
54647c6ae99SBarry Smith                                                          n3  n4  n5
54747c6ae99SBarry Smith                                                          n0  n1  n2
54847c6ae99SBarry Smith   */
54947c6ae99SBarry Smith 
55047c6ae99SBarry Smith   /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
55147c6ae99SBarry Smith   /* Assume Nodes are Internal to the Cube */
55247c6ae99SBarry Smith   n0  = rank - m*n - m - 1;
55347c6ae99SBarry Smith   n1  = rank - m*n - m;
55447c6ae99SBarry Smith   n2  = rank - m*n - m + 1;
55547c6ae99SBarry Smith   n3  = rank - m*n -1;
55647c6ae99SBarry Smith   n4  = rank - m*n;
55747c6ae99SBarry Smith   n5  = rank - m*n + 1;
55847c6ae99SBarry Smith   n6  = rank - m*n + m - 1;
55947c6ae99SBarry Smith   n7  = rank - m*n + m;
56047c6ae99SBarry Smith   n8  = rank - m*n + m + 1;
56147c6ae99SBarry Smith 
56247c6ae99SBarry Smith   n9  = rank - m - 1;
56347c6ae99SBarry Smith   n10 = rank - m;
56447c6ae99SBarry Smith   n11 = rank - m + 1;
56547c6ae99SBarry Smith   n12 = rank - 1;
56647c6ae99SBarry Smith   n14 = rank + 1;
56747c6ae99SBarry Smith   n15 = rank + m - 1;
56847c6ae99SBarry Smith   n16 = rank + m;
56947c6ae99SBarry Smith   n17 = rank + m + 1;
57047c6ae99SBarry Smith 
57147c6ae99SBarry Smith   n18 = rank + m*n - m - 1;
57247c6ae99SBarry Smith   n19 = rank + m*n - m;
57347c6ae99SBarry Smith   n20 = rank + m*n - m + 1;
57447c6ae99SBarry Smith   n21 = rank + m*n - 1;
57547c6ae99SBarry Smith   n22 = rank + m*n;
57647c6ae99SBarry Smith   n23 = rank + m*n + 1;
57747c6ae99SBarry Smith   n24 = rank + m*n + m - 1;
57847c6ae99SBarry Smith   n25 = rank + m*n + m;
57947c6ae99SBarry Smith   n26 = rank + m*n + m + 1;
58047c6ae99SBarry Smith 
58147c6ae99SBarry Smith   /* Assume Pieces are on Faces of Cube */
58247c6ae99SBarry Smith 
58347c6ae99SBarry Smith   if (xs == 0) { /* First assume not corner or edge */
58447c6ae99SBarry Smith     n0  = rank       -1 - (m*n);
58547c6ae99SBarry Smith     n3  = rank + m   -1 - (m*n);
58647c6ae99SBarry Smith     n6  = rank + 2*m -1 - (m*n);
58747c6ae99SBarry Smith     n9  = rank       -1;
58847c6ae99SBarry Smith     n12 = rank + m   -1;
58947c6ae99SBarry Smith     n15 = rank + 2*m -1;
59047c6ae99SBarry Smith     n18 = rank       -1 + (m*n);
59147c6ae99SBarry Smith     n21 = rank + m   -1 + (m*n);
59247c6ae99SBarry Smith     n24 = rank + 2*m -1 + (m*n);
59347c6ae99SBarry Smith    }
59447c6ae99SBarry Smith 
595ce00eea3SSatish Balay   if (xe == M) { /* First assume not corner or edge */
59647c6ae99SBarry Smith     n2  = rank -2*m +1 - (m*n);
59747c6ae99SBarry Smith     n5  = rank - m  +1 - (m*n);
59847c6ae99SBarry Smith     n8  = rank      +1 - (m*n);
59947c6ae99SBarry Smith     n11 = rank -2*m +1;
60047c6ae99SBarry Smith     n14 = rank - m  +1;
60147c6ae99SBarry Smith     n17 = rank      +1;
60247c6ae99SBarry Smith     n20 = rank -2*m +1 + (m*n);
60347c6ae99SBarry Smith     n23 = rank - m  +1 + (m*n);
60447c6ae99SBarry Smith     n26 = rank      +1 + (m*n);
60547c6ae99SBarry Smith   }
60647c6ae99SBarry Smith 
60747c6ae99SBarry Smith   if (ys==0) { /* First assume not corner or edge */
60847c6ae99SBarry Smith     n0  = rank + m * (n-1) -1 - (m*n);
60947c6ae99SBarry Smith     n1  = rank + m * (n-1)    - (m*n);
61047c6ae99SBarry Smith     n2  = rank + m * (n-1) +1 - (m*n);
61147c6ae99SBarry Smith     n9  = rank + m * (n-1) -1;
61247c6ae99SBarry Smith     n10 = rank + m * (n-1);
61347c6ae99SBarry Smith     n11 = rank + m * (n-1) +1;
61447c6ae99SBarry Smith     n18 = rank + m * (n-1) -1 + (m*n);
61547c6ae99SBarry Smith     n19 = rank + m * (n-1)    + (m*n);
61647c6ae99SBarry Smith     n20 = rank + m * (n-1) +1 + (m*n);
61747c6ae99SBarry Smith   }
61847c6ae99SBarry Smith 
61947c6ae99SBarry Smith   if (ye == N) { /* First assume not corner or edge */
62047c6ae99SBarry Smith     n6  = rank - m * (n-1) -1 - (m*n);
62147c6ae99SBarry Smith     n7  = rank - m * (n-1)    - (m*n);
62247c6ae99SBarry Smith     n8  = rank - m * (n-1) +1 - (m*n);
62347c6ae99SBarry Smith     n15 = rank - m * (n-1) -1;
62447c6ae99SBarry Smith     n16 = rank - m * (n-1);
62547c6ae99SBarry Smith     n17 = rank - m * (n-1) +1;
62647c6ae99SBarry Smith     n24 = rank - m * (n-1) -1 + (m*n);
62747c6ae99SBarry Smith     n25 = rank - m * (n-1)    + (m*n);
62847c6ae99SBarry Smith     n26 = rank - m * (n-1) +1 + (m*n);
62947c6ae99SBarry Smith   }
63047c6ae99SBarry Smith 
63147c6ae99SBarry Smith   if (zs == 0) { /* First assume not corner or edge */
63247c6ae99SBarry Smith     n0 = size - (m*n) + rank - m - 1;
63347c6ae99SBarry Smith     n1 = size - (m*n) + rank - m;
63447c6ae99SBarry Smith     n2 = size - (m*n) + rank - m + 1;
63547c6ae99SBarry Smith     n3 = size - (m*n) + rank - 1;
63647c6ae99SBarry Smith     n4 = size - (m*n) + rank;
63747c6ae99SBarry Smith     n5 = size - (m*n) + rank + 1;
63847c6ae99SBarry Smith     n6 = size - (m*n) + rank + m - 1;
63947c6ae99SBarry Smith     n7 = size - (m*n) + rank + m ;
64047c6ae99SBarry Smith     n8 = size - (m*n) + rank + m + 1;
64147c6ae99SBarry Smith   }
64247c6ae99SBarry Smith 
64347c6ae99SBarry Smith   if (ze == P) { /* First assume not corner or edge */
64447c6ae99SBarry Smith     n18 = (m*n) - (size-rank) - m - 1;
64547c6ae99SBarry Smith     n19 = (m*n) - (size-rank) - m;
64647c6ae99SBarry Smith     n20 = (m*n) - (size-rank) - m + 1;
64747c6ae99SBarry Smith     n21 = (m*n) - (size-rank) - 1;
64847c6ae99SBarry Smith     n22 = (m*n) - (size-rank);
64947c6ae99SBarry Smith     n23 = (m*n) - (size-rank) + 1;
65047c6ae99SBarry Smith     n24 = (m*n) - (size-rank) + m - 1;
65147c6ae99SBarry Smith     n25 = (m*n) - (size-rank) + m;
65247c6ae99SBarry Smith     n26 = (m*n) - (size-rank) + m + 1;
65347c6ae99SBarry Smith   }
65447c6ae99SBarry Smith 
65547c6ae99SBarry Smith   if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */
65647c6ae99SBarry Smith     n0 = size - m*n + rank + m-1 - m;
65747c6ae99SBarry Smith     n3 = size - m*n + rank + m-1;
65847c6ae99SBarry Smith     n6 = size - m*n + rank + m-1 + m;
65947c6ae99SBarry Smith   }
66047c6ae99SBarry Smith 
66147c6ae99SBarry Smith   if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */
66247c6ae99SBarry Smith     n18 = m*n - (size - rank) + m-1 - m;
66347c6ae99SBarry Smith     n21 = m*n - (size - rank) + m-1;
66447c6ae99SBarry Smith     n24 = m*n - (size - rank) + m-1 + m;
66547c6ae99SBarry Smith   }
66647c6ae99SBarry Smith 
66747c6ae99SBarry Smith   if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */
66847c6ae99SBarry Smith     n0  = rank + m*n -1 - m*n;
66947c6ae99SBarry Smith     n9  = rank + m*n -1;
67047c6ae99SBarry Smith     n18 = rank + m*n -1 + m*n;
67147c6ae99SBarry Smith   }
67247c6ae99SBarry Smith 
67347c6ae99SBarry Smith   if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */
67447c6ae99SBarry Smith     n6  = rank - m*(n-1) + m-1 - m*n;
67547c6ae99SBarry Smith     n15 = rank - m*(n-1) + m-1;
67647c6ae99SBarry Smith     n24 = rank - m*(n-1) + m-1 + m*n;
67747c6ae99SBarry Smith   }
67847c6ae99SBarry Smith 
679ce00eea3SSatish Balay   if ((xe==M) && (zs==0)) { /* Assume an edge, not corner */
68047c6ae99SBarry Smith     n2 = size - (m*n-rank) - (m-1) - m;
68147c6ae99SBarry Smith     n5 = size - (m*n-rank) - (m-1);
68247c6ae99SBarry Smith     n8 = size - (m*n-rank) - (m-1) + m;
68347c6ae99SBarry Smith   }
68447c6ae99SBarry Smith 
685ce00eea3SSatish Balay   if ((xe==M) && (ze==P)) { /* Assume an edge, not corner */
68647c6ae99SBarry Smith     n20 = m*n - (size - rank) - (m-1) - m;
68747c6ae99SBarry Smith     n23 = m*n - (size - rank) - (m-1);
68847c6ae99SBarry Smith     n26 = m*n - (size - rank) - (m-1) + m;
68947c6ae99SBarry Smith   }
69047c6ae99SBarry Smith 
691ce00eea3SSatish Balay   if ((xe==M) && (ys==0)) { /* Assume an edge, not corner */
69247c6ae99SBarry Smith     n2  = rank + m*(n-1) - (m-1) - m*n;
69347c6ae99SBarry Smith     n11 = rank + m*(n-1) - (m-1);
69447c6ae99SBarry Smith     n20 = rank + m*(n-1) - (m-1) + m*n;
69547c6ae99SBarry Smith   }
69647c6ae99SBarry Smith 
697ce00eea3SSatish Balay   if ((xe==M) && (ye==N)) { /* Assume an edge, not corner */
69847c6ae99SBarry Smith     n8  = rank - m*n +1 - m*n;
69947c6ae99SBarry Smith     n17 = rank - m*n +1;
70047c6ae99SBarry Smith     n26 = rank - m*n +1 + m*n;
70147c6ae99SBarry Smith   }
70247c6ae99SBarry Smith 
70347c6ae99SBarry Smith   if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */
70447c6ae99SBarry Smith     n0 = size - m + rank -1;
70547c6ae99SBarry Smith     n1 = size - m + rank;
70647c6ae99SBarry Smith     n2 = size - m + rank +1;
70747c6ae99SBarry Smith   }
70847c6ae99SBarry Smith 
70947c6ae99SBarry Smith   if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */
71047c6ae99SBarry Smith     n18 = m*n - (size - rank) + m*(n-1) -1;
71147c6ae99SBarry Smith     n19 = m*n - (size - rank) + m*(n-1);
71247c6ae99SBarry Smith     n20 = m*n - (size - rank) + m*(n-1) +1;
71347c6ae99SBarry Smith   }
71447c6ae99SBarry Smith 
71547c6ae99SBarry Smith   if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */
71647c6ae99SBarry Smith     n6 = size - (m*n-rank) - m * (n-1) -1;
71747c6ae99SBarry Smith     n7 = size - (m*n-rank) - m * (n-1);
71847c6ae99SBarry Smith     n8 = size - (m*n-rank) - m * (n-1) +1;
71947c6ae99SBarry Smith   }
72047c6ae99SBarry Smith 
72147c6ae99SBarry Smith   if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */
72247c6ae99SBarry Smith     n24 = rank - (size-m) -1;
72347c6ae99SBarry Smith     n25 = rank - (size-m);
72447c6ae99SBarry Smith     n26 = rank - (size-m) +1;
72547c6ae99SBarry Smith   }
72647c6ae99SBarry Smith 
72747c6ae99SBarry Smith   /* Check for Corners */
72847c6ae99SBarry Smith   if ((xs==0)   && (ys==0) && (zs==0)) { n0  = size -1;}
72947c6ae99SBarry Smith   if ((xs==0)   && (ys==0) && (ze==P)) { n18 = m*n-1;}
73047c6ae99SBarry Smith   if ((xs==0)   && (ye==N) && (zs==0)) { n6  = (size-1)-m*(n-1);}
73147c6ae99SBarry Smith   if ((xs==0)   && (ye==N) && (ze==P)) { n24 = m-1;}
732ce00eea3SSatish Balay   if ((xe==M) && (ys==0) && (zs==0)) { n2  = size-m;}
733ce00eea3SSatish Balay   if ((xe==M) && (ys==0) && (ze==P)) { n20 = m*n-m;}
734ce00eea3SSatish Balay   if ((xe==M) && (ye==N) && (zs==0)) { n8  = size-m*n;}
735ce00eea3SSatish Balay   if ((xe==M) && (ye==N) && (ze==P)) { n26 = 0;}
73647c6ae99SBarry Smith 
73747c6ae99SBarry Smith   /* Check for when not X,Y, and Z Periodic */
73847c6ae99SBarry Smith 
73947c6ae99SBarry Smith   /* If not X periodic */
7401321219cSEthan Coon   if (bx != DMDA_BOUNDARY_PERIODIC) {
74147c6ae99SBarry Smith     if (xs==0)   {n0  = n3  = n6  = n9  = n12 = n15 = n18 = n21 = n24 = -2;}
742ce00eea3SSatish Balay     if (xe==M) {n2  = n5  = n8  = n11 = n14 = n17 = n20 = n23 = n26 = -2;}
74347c6ae99SBarry Smith   }
74447c6ae99SBarry Smith 
74547c6ae99SBarry Smith   /* If not Y periodic */
7461321219cSEthan Coon   if (by != DMDA_BOUNDARY_PERIODIC) {
74747c6ae99SBarry Smith     if (ys==0)   {n0  = n1  = n2  = n9  = n10 = n11 = n18 = n19 = n20 = -2;}
74847c6ae99SBarry Smith     if (ye==N)   {n6  = n7  = n8  = n15 = n16 = n17 = n24 = n25 = n26 = -2;}
74947c6ae99SBarry Smith   }
75047c6ae99SBarry Smith 
75147c6ae99SBarry Smith   /* If not Z periodic */
7521321219cSEthan Coon   if (bz != DMDA_BOUNDARY_PERIODIC) {
75347c6ae99SBarry Smith     if (zs==0)   {n0  = n1  = n2  = n3  = n4  = n5  = n6  = n7  = n8  = -2;}
75447c6ae99SBarry Smith     if (ze==P)   {n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;}
75547c6ae99SBarry Smith   }
75647c6ae99SBarry Smith 
75747c6ae99SBarry Smith   ierr = PetscMalloc(27*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr);
75847c6ae99SBarry Smith   dd->neighbors[0] = n0;
75947c6ae99SBarry Smith   dd->neighbors[1] = n1;
76047c6ae99SBarry Smith   dd->neighbors[2] = n2;
76147c6ae99SBarry Smith   dd->neighbors[3] = n3;
76247c6ae99SBarry Smith   dd->neighbors[4] = n4;
76347c6ae99SBarry Smith   dd->neighbors[5] = n5;
76447c6ae99SBarry Smith   dd->neighbors[6] = n6;
76547c6ae99SBarry Smith   dd->neighbors[7] = n7;
76647c6ae99SBarry Smith   dd->neighbors[8] = n8;
76747c6ae99SBarry Smith   dd->neighbors[9] = n9;
76847c6ae99SBarry Smith   dd->neighbors[10] = n10;
76947c6ae99SBarry Smith   dd->neighbors[11] = n11;
77047c6ae99SBarry Smith   dd->neighbors[12] = n12;
77147c6ae99SBarry Smith   dd->neighbors[13] = rank;
77247c6ae99SBarry Smith   dd->neighbors[14] = n14;
77347c6ae99SBarry Smith   dd->neighbors[15] = n15;
77447c6ae99SBarry Smith   dd->neighbors[16] = n16;
77547c6ae99SBarry Smith   dd->neighbors[17] = n17;
77647c6ae99SBarry Smith   dd->neighbors[18] = n18;
77747c6ae99SBarry Smith   dd->neighbors[19] = n19;
77847c6ae99SBarry Smith   dd->neighbors[20] = n20;
77947c6ae99SBarry Smith   dd->neighbors[21] = n21;
78047c6ae99SBarry Smith   dd->neighbors[22] = n22;
78147c6ae99SBarry Smith   dd->neighbors[23] = n23;
78247c6ae99SBarry Smith   dd->neighbors[24] = n24;
78347c6ae99SBarry Smith   dd->neighbors[25] = n25;
78447c6ae99SBarry Smith   dd->neighbors[26] = n26;
78547c6ae99SBarry Smith 
78647c6ae99SBarry Smith   /* If star stencil then delete the corner neighbors */
78788661749SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR && o == 0) {
78847c6ae99SBarry Smith      /* save information about corner neighbors */
78947c6ae99SBarry Smith      sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7;
79047c6ae99SBarry Smith      sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18;
79147c6ae99SBarry Smith      sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25;
79247c6ae99SBarry Smith      sn26 = n26;
79347c6ae99SBarry Smith      n0  = n1  = n2  = n3  = n5  = n6  = n7  = n8  = n9  = n11 =
79447c6ae99SBarry Smith      n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1;
79547c6ae99SBarry Smith   }
79647c6ae99SBarry Smith 
79747c6ae99SBarry Smith 
79847c6ae99SBarry Smith   ierr = PetscMalloc((Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
79947c6ae99SBarry Smith   ierr = PetscLogObjectMemory(da,(Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt));CHKERRQ(ierr);
80047c6ae99SBarry Smith 
80147c6ae99SBarry Smith   nn = 0;
80247c6ae99SBarry Smith   /* Bottom Level */
80347c6ae99SBarry Smith   for (k=0; k<s_z; k++) {
80447c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
80547c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
806ce00eea3SSatish Balay         x_t = lx[n0 % m];
80747c6ae99SBarry Smith         y_t = ly[(n0 % (m*n))/m];
80847c6ae99SBarry Smith         z_t = lz[n0 / (m*n)];
80947c6ae99SBarry 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;
8106f951b95Secoon         if (twod && (s_t < 0)) {s_t = bases[n0] + x_t*y_t*z_t - (s_y-i)*x_t - s_x;} /* 2D case */
81147c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
81247c6ae99SBarry Smith       }
81347c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
81447c6ae99SBarry Smith         x_t = x;
81547c6ae99SBarry Smith         y_t = ly[(n1 % (m*n))/m];
81647c6ae99SBarry Smith         z_t = lz[n1 / (m*n)];
81747c6ae99SBarry 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;
8186f951b95Secoon         if (twod && (s_t < 0)) {s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t;} /* 2D case */
81947c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
82047c6ae99SBarry Smith       }
82147c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
822ce00eea3SSatish Balay         x_t = lx[n2 % m];
82347c6ae99SBarry Smith         y_t = ly[(n2 % (m*n))/m];
82447c6ae99SBarry Smith         z_t = lz[n2 / (m*n)];
82547c6ae99SBarry 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;
8266f951b95Secoon         if (twod && (s_t < 0)) {s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t;} /* 2D case */
82747c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
82847c6ae99SBarry Smith       }
82947c6ae99SBarry Smith     }
83047c6ae99SBarry Smith 
83147c6ae99SBarry Smith     for (i=0; i<y; i++) {
83247c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
833ce00eea3SSatish Balay         x_t = lx[n3 % m];
83447c6ae99SBarry Smith         y_t = y;
83547c6ae99SBarry Smith         z_t = lz[n3 / (m*n)];
83647c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8376f951b95Secoon         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 */
83847c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
83947c6ae99SBarry Smith       }
84047c6ae99SBarry Smith 
84147c6ae99SBarry Smith       if (n4 >= 0) { /* middle */
84247c6ae99SBarry Smith         x_t = x;
84347c6ae99SBarry Smith         y_t = y;
84447c6ae99SBarry Smith         z_t = lz[n4 / (m*n)];
84547c6ae99SBarry Smith         s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8466f951b95Secoon         if (twod && (s_t < 0)) {s_t = bases[n4] + i*x_t + x_t*y_t*z_t - x_t*y_t;} /* 2D case */
84747c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
848c2859e5eSBarry Smith       } else if (bz == DMDA_BOUNDARY_MIRROR) {
849c2859e5eSBarry Smith         for (j=0; j<x; j++) { idx[nn++] = 0;}
85047c6ae99SBarry Smith       }
85147c6ae99SBarry Smith 
85247c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
853ce00eea3SSatish Balay         x_t = lx[n5 % m];
85447c6ae99SBarry Smith         y_t = y;
85547c6ae99SBarry Smith         z_t = lz[n5 / (m*n)];
85647c6ae99SBarry Smith         s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8576f951b95Secoon         if (twod && (s_t < 0)) {s_t = bases[n5] + i*x_t + x_t*y_t*z_t - x_t*y_t;} /* 2D case */
85847c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
85947c6ae99SBarry Smith       }
86047c6ae99SBarry Smith     }
86147c6ae99SBarry Smith 
86247c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
86347c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
864ce00eea3SSatish Balay         x_t = lx[n6 % m];
86547c6ae99SBarry Smith         y_t = ly[(n6 % (m*n))/m];
86647c6ae99SBarry Smith         z_t = lz[n6 / (m*n)];
86747c6ae99SBarry Smith         s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8686f951b95Secoon         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 */
86947c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
87047c6ae99SBarry Smith       }
87147c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
87247c6ae99SBarry Smith         x_t = x;
87347c6ae99SBarry Smith         y_t = ly[(n7 % (m*n))/m];
87447c6ae99SBarry Smith         z_t = lz[n7 / (m*n)];
87547c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8766f951b95Secoon         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 */
87747c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
87847c6ae99SBarry Smith       }
87947c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
880ce00eea3SSatish Balay         x_t = lx[n8 % m];
88147c6ae99SBarry Smith         y_t = ly[(n8 % (m*n))/m];
88247c6ae99SBarry Smith         z_t = lz[n8 / (m*n)];
88347c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8846f951b95Secoon         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 */
88547c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
88647c6ae99SBarry Smith       }
88747c6ae99SBarry Smith     }
88847c6ae99SBarry Smith   }
88947c6ae99SBarry Smith 
89047c6ae99SBarry Smith   /* Middle Level */
89147c6ae99SBarry Smith   for (k=0; k<z; k++) {
89247c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
89347c6ae99SBarry Smith       if (n9 >= 0) { /* left below */
894ce00eea3SSatish Balay         x_t = lx[n9 % m];
89547c6ae99SBarry Smith         y_t = ly[(n9 % (m*n))/m];
89647c6ae99SBarry Smith         /* z_t = z; */
89747c6ae99SBarry Smith         s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
89847c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
89947c6ae99SBarry Smith       }
90047c6ae99SBarry Smith       if (n10 >= 0) { /* directly below */
90147c6ae99SBarry Smith         x_t = x;
90247c6ae99SBarry Smith         y_t = ly[(n10 % (m*n))/m];
90347c6ae99SBarry Smith         /* z_t = z; */
90447c6ae99SBarry Smith         s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
90547c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
906c2859e5eSBarry Smith       }  else if (by == DMDA_BOUNDARY_MIRROR) {
907c2859e5eSBarry Smith         for (j=0; j<x; j++) { idx[nn++] = 0;}
90847c6ae99SBarry Smith       }
90947c6ae99SBarry Smith       if (n11 >= 0) { /* right below */
910ce00eea3SSatish Balay         x_t = lx[n11 % m];
91147c6ae99SBarry Smith         y_t = ly[(n11 % (m*n))/m];
91247c6ae99SBarry Smith         /* z_t = z; */
91347c6ae99SBarry Smith         s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
91447c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
91547c6ae99SBarry Smith       }
91647c6ae99SBarry Smith     }
91747c6ae99SBarry Smith 
91847c6ae99SBarry Smith     for (i=0; i<y; i++) {
91947c6ae99SBarry Smith       if (n12 >= 0) { /* directly left */
920ce00eea3SSatish Balay         x_t = lx[n12 % m];
92147c6ae99SBarry Smith         y_t = y;
92247c6ae99SBarry Smith         /* z_t = z; */
92347c6ae99SBarry Smith         s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
92447c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
925c2859e5eSBarry Smith       }  else if (bx == DMDA_BOUNDARY_MIRROR) {
926c2859e5eSBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = 0;}
92747c6ae99SBarry Smith       }
92847c6ae99SBarry Smith 
92947c6ae99SBarry Smith       /* Interior */
93047c6ae99SBarry Smith       s_t = bases[rank] + i*x + k*x*y;
93147c6ae99SBarry Smith       for (j=0; j<x; j++) { idx[nn++] = s_t++;}
93247c6ae99SBarry Smith 
93347c6ae99SBarry Smith       if (n14 >= 0) { /* directly right */
934ce00eea3SSatish Balay         x_t = lx[n14 % m];
93547c6ae99SBarry Smith         y_t = y;
93647c6ae99SBarry Smith         /* z_t = z; */
93747c6ae99SBarry Smith         s_t = bases[n14] + i*x_t + k*x_t*y_t;
93847c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
939c2859e5eSBarry Smith       } else if (bx == DMDA_BOUNDARY_MIRROR) {
940c2859e5eSBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = 0;}
94147c6ae99SBarry Smith       }
94247c6ae99SBarry Smith     }
94347c6ae99SBarry Smith 
94447c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
94547c6ae99SBarry Smith       if (n15 >= 0) { /* left above */
946ce00eea3SSatish Balay         x_t = lx[n15 % m];
94747c6ae99SBarry Smith         y_t = ly[(n15 % (m*n))/m];
94847c6ae99SBarry Smith         /* z_t = z; */
94947c6ae99SBarry Smith         s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
95047c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
95147c6ae99SBarry Smith       }
95247c6ae99SBarry Smith       if (n16 >= 0) { /* directly above */
95347c6ae99SBarry Smith         x_t = x;
95447c6ae99SBarry Smith         y_t = ly[(n16 % (m*n))/m];
95547c6ae99SBarry Smith         /* z_t = z; */
95647c6ae99SBarry Smith         s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
95747c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
958c2859e5eSBarry Smith       } else if (by == DMDA_BOUNDARY_MIRROR) {
959c2859e5eSBarry Smith         for (j=0; j<x; j++) { idx[nn++] = 0;}
96047c6ae99SBarry Smith       }
96147c6ae99SBarry Smith       if (n17 >= 0) { /* right above */
962ce00eea3SSatish Balay         x_t = lx[n17 % m];
96347c6ae99SBarry Smith         y_t = ly[(n17 % (m*n))/m];
96447c6ae99SBarry Smith         /* z_t = z; */
96547c6ae99SBarry Smith         s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
96647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
96747c6ae99SBarry Smith       }
96847c6ae99SBarry Smith     }
96947c6ae99SBarry Smith   }
97047c6ae99SBarry Smith 
97147c6ae99SBarry Smith   /* Upper Level */
97247c6ae99SBarry Smith   for (k=0; k<s_z; k++) {
97347c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
97447c6ae99SBarry Smith       if (n18 >= 0) { /* left below */
975ce00eea3SSatish Balay         x_t = lx[n18 % m];
97647c6ae99SBarry Smith         y_t = ly[(n18 % (m*n))/m];
97747c6ae99SBarry Smith         /* z_t = lz[n18 / (m*n)]; */
97847c6ae99SBarry Smith         s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
9796f951b95Secoon         if (twod && (s_t >= M*N*P)) {s_t = bases[n18] - (s_y-i)*x_t -s_x + x_t*y_t;} /* 2d case */
98047c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
98147c6ae99SBarry Smith       }
98247c6ae99SBarry Smith       if (n19 >= 0) { /* directly below */
98347c6ae99SBarry Smith         x_t = x;
98447c6ae99SBarry Smith         y_t = ly[(n19 % (m*n))/m];
98547c6ae99SBarry Smith         /* z_t = lz[n19 / (m*n)]; */
98647c6ae99SBarry Smith         s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
9876f951b95Secoon         if (twod && (s_t >= M*N*P)) {s_t = bases[n19] - (s_y+1-i)*x_t + x_t*y_t;} /* 2d case */
98847c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
98947c6ae99SBarry Smith       }
99047c6ae99SBarry Smith       if (n20 >= 0) { /* right below */
991ce00eea3SSatish Balay         x_t = lx[n20 % m];
99247c6ae99SBarry Smith         y_t = ly[(n20 % (m*n))/m];
99347c6ae99SBarry Smith         /* z_t = lz[n20 / (m*n)]; */
99447c6ae99SBarry Smith         s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
9956f951b95Secoon         if (twod && (s_t >= M*N*P)) {s_t = bases[n20] - (s_y+1-i)*x_t + x_t*y_t;} /* 2d case */
99647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
99747c6ae99SBarry Smith       }
99847c6ae99SBarry Smith     }
99947c6ae99SBarry Smith 
100047c6ae99SBarry Smith     for (i=0; i<y; i++) {
100147c6ae99SBarry Smith       if (n21 >= 0) { /* directly left */
1002ce00eea3SSatish Balay         x_t = lx[n21 % m];
100347c6ae99SBarry Smith         y_t = y;
100447c6ae99SBarry Smith         /* z_t = lz[n21 / (m*n)]; */
100547c6ae99SBarry Smith         s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
10066f951b95Secoon         if (twod && (s_t >= M*N*P)) {s_t = bases[n21] + (i+1)*x_t - s_x;}  /* 2d case */
100747c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
100847c6ae99SBarry Smith       }
100947c6ae99SBarry Smith 
101047c6ae99SBarry Smith       if (n22 >= 0) { /* middle */
101147c6ae99SBarry Smith         x_t = x;
101247c6ae99SBarry Smith         y_t = y;
101347c6ae99SBarry Smith         /* z_t = lz[n22 / (m*n)]; */
101447c6ae99SBarry Smith         s_t = bases[n22] + i*x_t + k*x_t*y_t;
10156f951b95Secoon         if (twod && (s_t >= M*N*P)) {s_t = bases[n22] + i*x_t;} /* 2d case */
101647c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1017c2859e5eSBarry Smith       } else if (bz == DMDA_BOUNDARY_MIRROR) {
1018c2859e5eSBarry Smith         for (j=0; j<x; j++) { idx[nn++] = 0;}
101947c6ae99SBarry Smith       }
102047c6ae99SBarry Smith 
102147c6ae99SBarry Smith       if (n23 >= 0) { /* directly right */
1022ce00eea3SSatish Balay         x_t = lx[n23 % m];
102347c6ae99SBarry Smith         y_t = y;
102447c6ae99SBarry Smith         /* z_t = lz[n23 / (m*n)]; */
102547c6ae99SBarry Smith         s_t = bases[n23] + i*x_t + k*x_t*y_t;
10266f951b95Secoon         if (twod && (s_t >= M*N*P)) {s_t = bases[n23] + i*x_t;} /* 2d case */
102747c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
102847c6ae99SBarry Smith       }
102947c6ae99SBarry Smith     }
103047c6ae99SBarry Smith 
103147c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
103247c6ae99SBarry Smith       if (n24 >= 0) { /* left above */
1033ce00eea3SSatish Balay         x_t = lx[n24 % m];
103447c6ae99SBarry Smith         y_t = ly[(n24 % (m*n))/m];
103547c6ae99SBarry Smith         /* z_t = lz[n24 / (m*n)]; */
103647c6ae99SBarry Smith         s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
10376f951b95Secoon         if (twod && (s_t >= M*N*P)) {s_t = bases[n24] + i*x_t - s_x;} /* 2d case */
103847c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
103947c6ae99SBarry Smith       }
104047c6ae99SBarry Smith       if (n25 >= 0) { /* directly above */
104147c6ae99SBarry Smith         x_t = x;
104247c6ae99SBarry Smith         y_t = ly[(n25 % (m*n))/m];
104347c6ae99SBarry Smith         /* z_t = lz[n25 / (m*n)]; */
104447c6ae99SBarry Smith         s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
10456f951b95Secoon         if (twod && (s_t >= M*N*P)) {s_t = bases[n25] + (i-1)*x_t;} /* 2d case */
104647c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
104747c6ae99SBarry Smith       }
104847c6ae99SBarry Smith       if (n26 >= 0) { /* right above */
1049ce00eea3SSatish Balay         x_t = lx[n26 % m];
105047c6ae99SBarry Smith         y_t = ly[(n26 % (m*n))/m];
105147c6ae99SBarry Smith         /* z_t = lz[n26 / (m*n)]; */
105247c6ae99SBarry Smith         s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
10536f951b95Secoon         if (twod && (s_t >= M*N*P)) {s_t = bases[n26] + (i-1)*x_t;} /* 2d case */
105447c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
105547c6ae99SBarry Smith       }
105647c6ae99SBarry Smith     }
105747c6ae99SBarry Smith   }
1058ce00eea3SSatish Balay 
1059ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
106047c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
106147c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr);
1062fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
1063fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
106447c6ae99SBarry Smith 
106588661749SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR && o == 0) {
106647c6ae99SBarry Smith     n0  = sn0;  n1  = sn1;  n2  = sn2;  n3  = sn3;  n5  = sn5;  n6  = sn6; n7 = sn7;
106747c6ae99SBarry Smith     n8  = sn8;  n9  = sn9;  n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18;
106847c6ae99SBarry Smith     n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25;
106947c6ae99SBarry Smith     n26 = sn26;
1070ce00eea3SSatish Balay   }
107147c6ae99SBarry Smith 
107288661749SPeter Brune   if (((stencil_type == DMDA_STENCIL_STAR) ||
10731321219cSEthan Coon       (bx != DMDA_BOUNDARY_PERIODIC && bx) ||
10741321219cSEthan Coon       (by != DMDA_BOUNDARY_PERIODIC && by) ||
107588661749SPeter Brune        (bz != DMDA_BOUNDARY_PERIODIC && bz)) && o == 0) {
1076ce00eea3SSatish Balay     /*
1077ce00eea3SSatish Balay         Recompute the local to global mappings, this time keeping the
1078ce00eea3SSatish Balay       information about the cross corner processor numbers.
1079ce00eea3SSatish Balay     */
108047c6ae99SBarry Smith     nn = 0;
108147c6ae99SBarry Smith     /* Bottom Level */
108247c6ae99SBarry Smith     for (k=0; k<s_z; k++) {
108347c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
108447c6ae99SBarry Smith         if (n0 >= 0) { /* left below */
1085ce00eea3SSatish Balay           x_t = lx[n0 % m];
108647c6ae99SBarry Smith           y_t = ly[(n0 % (m*n))/m];
108747c6ae99SBarry Smith           z_t = lz[n0 / (m*n)];
108847c6ae99SBarry 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;
108947c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1090ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Ys-ys < 0 && Zs-zs < 0) {
1091ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
109247c6ae99SBarry Smith         }
109347c6ae99SBarry Smith         if (n1 >= 0) { /* directly below */
109447c6ae99SBarry Smith           x_t = x;
109547c6ae99SBarry Smith           y_t = ly[(n1 % (m*n))/m];
109647c6ae99SBarry Smith           z_t = lz[n1 / (m*n)];
109747c6ae99SBarry 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;
109847c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1099ce00eea3SSatish Balay         } else if (Ys-ys < 0 && Zs-zs < 0) {
1100ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
110147c6ae99SBarry Smith         }
110247c6ae99SBarry Smith         if (n2 >= 0) { /* right below */
1103ce00eea3SSatish Balay           x_t = lx[n2 % m];
110447c6ae99SBarry Smith           y_t = ly[(n2 % (m*n))/m];
110547c6ae99SBarry Smith           z_t = lz[n2 / (m*n)];
110647c6ae99SBarry 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;
110747c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1108ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Ys-ys < 0 && Zs-zs < 0) {
1109ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
111047c6ae99SBarry Smith         }
111147c6ae99SBarry Smith       }
111247c6ae99SBarry Smith 
111347c6ae99SBarry Smith       for (i=0; i<y; i++) {
111447c6ae99SBarry Smith         if (n3 >= 0) { /* directly left */
1115ce00eea3SSatish Balay           x_t = lx[n3 % m];
111647c6ae99SBarry Smith           y_t = y;
111747c6ae99SBarry Smith           z_t = lz[n3 / (m*n)];
111847c6ae99SBarry Smith           s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
111947c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1120ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Zs-zs < 0) {
1121ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
112247c6ae99SBarry Smith         }
112347c6ae99SBarry Smith 
112447c6ae99SBarry Smith         if (n4 >= 0) { /* middle */
112547c6ae99SBarry Smith           x_t = x;
112647c6ae99SBarry Smith           y_t = y;
112747c6ae99SBarry Smith           z_t = lz[n4 / (m*n)];
112847c6ae99SBarry Smith           s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
112947c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1130ce00eea3SSatish Balay         } else if (Zs-zs < 0) {
1131c2859e5eSBarry Smith           if (bz == DMDA_BOUNDARY_MIRROR) {
1132c2859e5eSBarry Smith             for (j=0; j<x; j++) { idx[nn++] = 0;}
1133c2859e5eSBarry Smith           } else {
1134ce00eea3SSatish Balay             for (j=0; j<x; j++) { idx[nn++] = -1;}
113547c6ae99SBarry Smith           }
1136c2859e5eSBarry Smith         }
113747c6ae99SBarry Smith 
113847c6ae99SBarry Smith         if (n5 >= 0) { /* directly right */
1139ce00eea3SSatish Balay           x_t = lx[n5 % m];
114047c6ae99SBarry Smith           y_t = y;
114147c6ae99SBarry Smith           z_t = lz[n5 / (m*n)];
114247c6ae99SBarry Smith           s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
114347c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1144ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Zs-zs < 0) {
1145ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
114647c6ae99SBarry Smith         }
114747c6ae99SBarry Smith       }
114847c6ae99SBarry Smith 
114947c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
115047c6ae99SBarry Smith         if (n6 >= 0) { /* left above */
1151ce00eea3SSatish Balay           x_t = lx[n6 % m];
115247c6ae99SBarry Smith           y_t = ly[(n6 % (m*n))/m];
115347c6ae99SBarry Smith           z_t = lz[n6 / (m*n)];
115447c6ae99SBarry Smith           s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
115547c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1156ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ye-Ye < 0 && Zs-zs < 0) {
1157ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
115847c6ae99SBarry Smith         }
115947c6ae99SBarry Smith         if (n7 >= 0) { /* directly above */
116047c6ae99SBarry Smith           x_t = x;
116147c6ae99SBarry Smith           y_t = ly[(n7 % (m*n))/m];
116247c6ae99SBarry Smith           z_t = lz[n7 / (m*n)];
116347c6ae99SBarry Smith           s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
116447c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1165ce00eea3SSatish Balay         } else if (ye-Ye < 0 && Zs-zs < 0) {
1166ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
116747c6ae99SBarry Smith         }
116847c6ae99SBarry Smith         if (n8 >= 0) { /* right above */
1169ce00eea3SSatish Balay           x_t = lx[n8 % m];
117047c6ae99SBarry Smith           y_t = ly[(n8 % (m*n))/m];
117147c6ae99SBarry Smith           z_t = lz[n8 / (m*n)];
117247c6ae99SBarry Smith           s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
117347c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1174ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ye-Ye < 0 && Zs-zs < 0) {
1175ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
117647c6ae99SBarry Smith         }
117747c6ae99SBarry Smith       }
117847c6ae99SBarry Smith     }
117947c6ae99SBarry Smith 
118047c6ae99SBarry Smith     /* Middle Level */
118147c6ae99SBarry Smith     for (k=0; k<z; k++) {
118247c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
118347c6ae99SBarry Smith         if (n9 >= 0) { /* left below */
1184ce00eea3SSatish Balay           x_t = lx[n9 % m];
118547c6ae99SBarry Smith           y_t = ly[(n9 % (m*n))/m];
118647c6ae99SBarry Smith           /* z_t = z; */
118747c6ae99SBarry Smith           s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
118847c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1189ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Ys-ys < 0) {
1190ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
119147c6ae99SBarry Smith         }
119247c6ae99SBarry Smith         if (n10 >= 0) { /* directly below */
119347c6ae99SBarry Smith           x_t = x;
119447c6ae99SBarry Smith           y_t = ly[(n10 % (m*n))/m];
119547c6ae99SBarry Smith           /* z_t = z; */
119647c6ae99SBarry Smith           s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
119747c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1198ce00eea3SSatish Balay         } else if (Ys-ys < 0) {
1199c2859e5eSBarry Smith           if (by == DMDA_BOUNDARY_MIRROR) {
1200ce00eea3SSatish Balay             for (j=0; j<x; j++) { idx[nn++] = -1;}
1201c2859e5eSBarry Smith           } else {
1202c2859e5eSBarry Smith             for (j=0; j<x; j++) { idx[nn++] = -1;}
1203c2859e5eSBarry Smith           }
120447c6ae99SBarry Smith         }
120547c6ae99SBarry Smith         if (n11 >= 0) { /* right below */
1206ce00eea3SSatish Balay           x_t = lx[n11 % m];
120747c6ae99SBarry Smith           y_t = ly[(n11 % (m*n))/m];
120847c6ae99SBarry Smith           /* z_t = z; */
120947c6ae99SBarry Smith           s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
121047c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1211ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Ys-ys < 0) {
1212ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
121347c6ae99SBarry Smith         }
121447c6ae99SBarry Smith       }
121547c6ae99SBarry Smith 
121647c6ae99SBarry Smith       for (i=0; i<y; i++) {
121747c6ae99SBarry Smith         if (n12 >= 0) { /* directly left */
1218ce00eea3SSatish Balay           x_t = lx[n12 % m];
121947c6ae99SBarry Smith           y_t = y;
122047c6ae99SBarry Smith           /* z_t = z; */
122147c6ae99SBarry Smith           s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
122247c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1223ce00eea3SSatish Balay         } else if (Xs-xs < 0) {
1224c2859e5eSBarry Smith           if (bx == DMDA_BOUNDARY_MIRROR) {
1225c2859e5eSBarry Smith             for (j=0; j<s_x; j++) { idx[nn++] = 0;}
1226c2859e5eSBarry Smith           } else {
1227ce00eea3SSatish Balay             for (j=0; j<s_x; j++) { idx[nn++] = -1;}
122847c6ae99SBarry Smith           }
1229c2859e5eSBarry Smith         }
123047c6ae99SBarry Smith 
123147c6ae99SBarry Smith         /* Interior */
123247c6ae99SBarry Smith         s_t = bases[rank] + i*x + k*x*y;
123347c6ae99SBarry Smith         for (j=0; j<x; j++) { idx[nn++] = s_t++;}
123447c6ae99SBarry Smith 
123547c6ae99SBarry Smith         if (n14 >= 0) { /* directly right */
1236ce00eea3SSatish Balay           x_t = lx[n14 % m];
123747c6ae99SBarry Smith           y_t = y;
123847c6ae99SBarry Smith           /* z_t = z; */
123947c6ae99SBarry Smith           s_t = bases[n14] + i*x_t + k*x_t*y_t;
124047c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1241ce00eea3SSatish Balay         } else if (xe-Xe < 0) {
1242c2859e5eSBarry Smith           if (bx == DMDA_BOUNDARY_MIRROR) {
1243c2859e5eSBarry Smith             for (j=0; j<s_x; j++) { idx[nn++] = 0;}
1244c2859e5eSBarry Smith           } else {
1245ce00eea3SSatish Balay             for (j=0; j<s_x; j++) { idx[nn++] = -1;}
124647c6ae99SBarry Smith           }
124747c6ae99SBarry Smith         }
1248c2859e5eSBarry Smith       }
124947c6ae99SBarry Smith 
125047c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
125147c6ae99SBarry Smith         if (n15 >= 0) { /* left above */
1252ce00eea3SSatish Balay           x_t = lx[n15 % m];
125347c6ae99SBarry Smith           y_t = ly[(n15 % (m*n))/m];
125447c6ae99SBarry Smith           /* z_t = z; */
125547c6ae99SBarry Smith           s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
125647c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1257ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ye-Ye < 0) {
1258ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
125947c6ae99SBarry Smith         }
126047c6ae99SBarry Smith         if (n16 >= 0) { /* directly above */
126147c6ae99SBarry Smith           x_t = x;
126247c6ae99SBarry Smith           y_t = ly[(n16 % (m*n))/m];
126347c6ae99SBarry Smith           /* z_t = z; */
126447c6ae99SBarry Smith           s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
126547c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1266ce00eea3SSatish Balay         } else if (ye-Ye < 0) {
1267c2859e5eSBarry Smith           if (by == DMDA_BOUNDARY_MIRROR) {
1268c2859e5eSBarry Smith             for (j=0; j<x; j++) { idx[nn++] = 0;}
1269c2859e5eSBarry Smith           } else {
1270ce00eea3SSatish Balay             for (j=0; j<x; j++) { idx[nn++] = -1;}
127147c6ae99SBarry Smith           }
1272c2859e5eSBarry Smith         }
127347c6ae99SBarry Smith         if (n17 >= 0) { /* right above */
1274ce00eea3SSatish Balay           x_t = lx[n17 % m];
127547c6ae99SBarry Smith           y_t = ly[(n17 % (m*n))/m];
127647c6ae99SBarry Smith           /* z_t = z; */
127747c6ae99SBarry Smith           s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
127847c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1279ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ye-Ye < 0) {
1280ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
128147c6ae99SBarry Smith         }
128247c6ae99SBarry Smith       }
128347c6ae99SBarry Smith     }
128447c6ae99SBarry Smith 
128547c6ae99SBarry Smith     /* Upper Level */
128647c6ae99SBarry Smith     for (k=0; k<s_z; k++) {
128747c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
128847c6ae99SBarry Smith         if (n18 >= 0) { /* left below */
1289ce00eea3SSatish Balay           x_t = lx[n18 % m];
129047c6ae99SBarry Smith           y_t = ly[(n18 % (m*n))/m];
129147c6ae99SBarry Smith           /* z_t = lz[n18 / (m*n)]; */
129247c6ae99SBarry Smith           s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
129347c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1294ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Ys-ys < 0 && ze-Ze < 0) {
1295ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
129647c6ae99SBarry Smith         }
129747c6ae99SBarry Smith         if (n19 >= 0) { /* directly below */
129847c6ae99SBarry Smith           x_t = x;
129947c6ae99SBarry Smith           y_t = ly[(n19 % (m*n))/m];
130047c6ae99SBarry Smith           /* z_t = lz[n19 / (m*n)]; */
130147c6ae99SBarry Smith           s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
130247c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1303ce00eea3SSatish Balay         } else if (Ys-ys < 0 && ze-Ze < 0) {
1304ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
130547c6ae99SBarry Smith         }
130647c6ae99SBarry Smith         if (n20 >= 0) { /* right below */
1307ce00eea3SSatish Balay           x_t = lx[n20 % m];
130847c6ae99SBarry Smith           y_t = ly[(n20 % (m*n))/m];
130947c6ae99SBarry Smith           /* z_t = lz[n20 / (m*n)]; */
131047c6ae99SBarry Smith           s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
131147c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1312ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Ys-ys < 0 && ze-Ze < 0) {
1313ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
131447c6ae99SBarry Smith         }
131547c6ae99SBarry Smith       }
131647c6ae99SBarry Smith 
131747c6ae99SBarry Smith       for (i=0; i<y; i++) {
131847c6ae99SBarry Smith         if (n21 >= 0) { /* directly left */
1319ce00eea3SSatish Balay           x_t = lx[n21 % m];
132047c6ae99SBarry Smith           y_t = y;
132147c6ae99SBarry Smith           /* z_t = lz[n21 / (m*n)]; */
132247c6ae99SBarry Smith           s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
132347c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1324ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ze-Ze < 0) {
1325ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
132647c6ae99SBarry Smith         }
132747c6ae99SBarry Smith 
132847c6ae99SBarry Smith         if (n22 >= 0) { /* middle */
132947c6ae99SBarry Smith           x_t = x;
133047c6ae99SBarry Smith           y_t = y;
133147c6ae99SBarry Smith           /* z_t = lz[n22 / (m*n)]; */
133247c6ae99SBarry Smith           s_t = bases[n22] + i*x_t + k*x_t*y_t;
133347c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1334ce00eea3SSatish Balay         } else if (ze-Ze < 0) {
1335c2859e5eSBarry Smith           if (bz == DMDA_BOUNDARY_MIRROR) {
1336c2859e5eSBarry Smith             for (j=0; j<x; j++) { idx[nn++] = 0;}
1337c2859e5eSBarry Smith           } else {
1338ce00eea3SSatish Balay             for (j=0; j<x; j++) { idx[nn++] = -1;}
133947c6ae99SBarry Smith           }
1340c2859e5eSBarry Smith         }
134147c6ae99SBarry Smith 
134247c6ae99SBarry Smith         if (n23 >= 0) { /* directly right */
1343ce00eea3SSatish Balay           x_t = lx[n23 % m];
134447c6ae99SBarry Smith           y_t = y;
134547c6ae99SBarry Smith           /* z_t = lz[n23 / (m*n)]; */
134647c6ae99SBarry Smith           s_t = bases[n23] + i*x_t + k*x_t*y_t;
134747c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1348ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ze-Ze < 0) {
1349ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
135047c6ae99SBarry Smith         }
135147c6ae99SBarry Smith       }
135247c6ae99SBarry Smith 
135347c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
135447c6ae99SBarry Smith         if (n24 >= 0) { /* left above */
1355ce00eea3SSatish Balay           x_t = lx[n24 % m];
135647c6ae99SBarry Smith           y_t = ly[(n24 % (m*n))/m];
135747c6ae99SBarry Smith           /* z_t = lz[n24 / (m*n)]; */
135847c6ae99SBarry Smith           s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
135947c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1360ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ye-Ye < 0 && ze-Ze < 0) {
1361ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
136247c6ae99SBarry Smith         }
136347c6ae99SBarry Smith         if (n25 >= 0) { /* directly above */
136447c6ae99SBarry Smith           x_t = x;
136547c6ae99SBarry Smith           y_t = ly[(n25 % (m*n))/m];
136647c6ae99SBarry Smith           /* z_t = lz[n25 / (m*n)]; */
136747c6ae99SBarry Smith           s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
136847c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1369ce00eea3SSatish Balay         } else if (ye-Ye < 0 && ze-Ze < 0) {
1370ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
137147c6ae99SBarry Smith         }
137247c6ae99SBarry Smith         if (n26 >= 0) { /* right above */
1373ce00eea3SSatish Balay           x_t = lx[n26 % m];
137447c6ae99SBarry Smith           y_t = ly[(n26 % (m*n))/m];
137547c6ae99SBarry Smith           /* z_t = lz[n26 / (m*n)]; */
137647c6ae99SBarry Smith           s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
137747c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1378ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ye-Ye < 0 && ze-Ze < 0) {
1379ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
138047c6ae99SBarry Smith         }
138147c6ae99SBarry Smith       }
138247c6ae99SBarry Smith     }
138347c6ae99SBarry Smith   }
138447c6ae99SBarry Smith   /*
138547c6ae99SBarry Smith      Set the local to global ordering in the global vector, this allows use
138647c6ae99SBarry Smith      of VecSetValuesLocal().
138747c6ae99SBarry Smith   */
1388ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_OWN_POINTER,&ltogis);CHKERRQ(ierr);
1389ce00eea3SSatish Balay   ierr = PetscMalloc(nn*dof*sizeof(PetscInt),&idx_cpy);CHKERRQ(ierr);
1390db87c5ecSEthan Coon   ierr = PetscLogObjectMemory(da,nn*dof*sizeof(PetscInt));CHKERRQ(ierr);
1391*3bf036e2SBarry Smith   ierr = ISGetIndices(ltogis, &idx_full);CHKERRQ(ierr);
1392ce00eea3SSatish Balay   ierr = PetscMemcpy(idx_cpy,idx_full,nn*dof*sizeof(PetscInt));CHKERRQ(ierr);
1393*3bf036e2SBarry Smith   ierr = ISRestoreIndices(ltogis, &idx_full);CHKERRQ(ierr);
1394ce00eea3SSatish Balay   ierr = ISLocalToGlobalMappingCreateIS(ltogis,&da->ltogmap);CHKERRQ(ierr);
1395ce00eea3SSatish Balay   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
1396fcfd50ebSBarry Smith   ierr = ISDestroy(&ltogis);CHKERRQ(ierr);
13971411c6eeSJed Brown   ierr = ISLocalToGlobalMappingBlock(da->ltogmap,dd->w,&da->ltogmapb);CHKERRQ(ierr);
13981411c6eeSJed Brown   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
139947c6ae99SBarry Smith 
1400ce00eea3SSatish Balay   ierr = PetscFree2(bases,ldims);CHKERRQ(ierr);
1401ce00eea3SSatish Balay   dd->m  = m;  dd->n  = n;  dd->p  = p;
1402ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
1403ce00eea3SSatish Balay   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = zs; dd->ze = ze;
1404ce00eea3SSatish Balay   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = Zs; dd->Ze = Ze;
1405ce00eea3SSatish Balay 
1406fcfd50ebSBarry Smith   ierr = VecDestroy(&local);CHKERRQ(ierr);
1407fcfd50ebSBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
1408ce00eea3SSatish Balay 
1409ce00eea3SSatish Balay   dd->gtol      = gtol;
1410ce00eea3SSatish Balay   dd->ltog      = ltog;
1411ce00eea3SSatish Balay   dd->idx       = idx_cpy;
1412ce00eea3SSatish Balay   dd->Nl        = nn*dof;
1413ce00eea3SSatish Balay   dd->base      = base;
1414ce00eea3SSatish Balay   da->ops->view = DMView_DA_3d;
141547c6ae99SBarry Smith   dd->ltol = PETSC_NULL;
141647c6ae99SBarry Smith   dd->ao   = PETSC_NULL;
1417ce00eea3SSatish Balay 
141847c6ae99SBarry Smith   PetscFunctionReturn(0);
141947c6ae99SBarry Smith }
1420564755cdSBarry Smith 
142147c6ae99SBarry Smith 
142247c6ae99SBarry Smith #undef __FUNCT__
1423aa219208SBarry Smith #define __FUNCT__ "DMDACreate3d"
142447c6ae99SBarry Smith /*@C
1425aa219208SBarry Smith    DMDACreate3d - Creates an object that will manage the communication of three-dimensional
142647c6ae99SBarry Smith    regular array data that is distributed across some processors.
142747c6ae99SBarry Smith 
142847c6ae99SBarry Smith    Collective on MPI_Comm
142947c6ae99SBarry Smith 
143047c6ae99SBarry Smith    Input Parameters:
143147c6ae99SBarry Smith +  comm - MPI communicator
14321321219cSEthan Coon .  bx,by,bz - type of ghost nodes the array have.
14331321219cSEthan Coon          Use one of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC.
1434aa219208SBarry Smith .  stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX)
143547c6ae99SBarry 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
143647c6ae99SBarry Smith             from the command line with -da_grid_x <M> -da_grid_y <N> -da_grid_z <P>)
143747c6ae99SBarry Smith .  m,n,p - corresponding number of processors in each dimension
143847c6ae99SBarry Smith            (or PETSC_DECIDE to have calculated)
143947c6ae99SBarry Smith .  dof - number of degrees of freedom per node
144010d7c030SMatthew G Knepley .  s - stencil width
144110d7c030SMatthew G Knepley -  lx, ly, lz - arrays containing the number of nodes in each cell along
144247c6ae99SBarry Smith           the x, y, and z coordinates, or PETSC_NULL. If non-null, these
144347c6ae99SBarry Smith           must be of length as m,n,p and the corresponding
144447c6ae99SBarry Smith           m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of
144547c6ae99SBarry Smith           the ly[] must N, sum of the lz[] must be P
144647c6ae99SBarry Smith 
144747c6ae99SBarry Smith    Output Parameter:
144847c6ae99SBarry Smith .  da - the resulting distributed array object
144947c6ae99SBarry Smith 
145047c6ae99SBarry Smith    Options Database Key:
1451706a11cbSBarry Smith +  -dm_view - Calls DMView() at the conclusion of DMDACreate3d()
145247c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
145347c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
145447c6ae99SBarry Smith .  -da_grid_z <nz> - number of grid points in z direction, if P < 0
145547c6ae99SBarry Smith .  -da_processors_x <MX> - number of processors in x direction
145647c6ae99SBarry Smith .  -da_processors_y <MY> - number of processors in y direction
145747c6ae99SBarry Smith .  -da_processors_z <MZ> - number of processors in z direction
1458e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement ratio in x direction
1459e0f5d30fSBarry Smith .  -da_refine_y <ry> - refinement ratio in y direction
1460e0f5d30fSBarry Smith .  -da_refine_z <rz>- refinement ratio in z directio
1461e0f5d30fSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating it, , if M, N, or P < 0
146247c6ae99SBarry Smith 
146347c6ae99SBarry Smith    Level: beginner
146447c6ae99SBarry Smith 
146547c6ae99SBarry Smith    Notes:
1466aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
1467aa219208SBarry Smith    standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
146847c6ae99SBarry Smith    the standard 27-pt stencil.
146947c6ae99SBarry Smith 
1470aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
1471564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
1472564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
147347c6ae99SBarry Smith 
147447c6ae99SBarry Smith .keywords: distributed array, create, three-dimensional
147547c6ae99SBarry Smith 
1476aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate2d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
1477aa219208SBarry Smith           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDASetRefinementFactor(),
1478d461ba97SBarry Smith           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
147947c6ae99SBarry Smith 
148047c6ae99SBarry Smith @*/
14811321219cSEthan Coon PetscErrorCode  DMDACreate3d(MPI_Comm comm,DMDABoundaryType bx,DMDABoundaryType by,DMDABoundaryType bz,DMDAStencilType stencil_type,PetscInt M,
14829a42bb27SBarry 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)
148347c6ae99SBarry Smith {
148447c6ae99SBarry Smith   PetscErrorCode ierr;
148547c6ae99SBarry Smith 
148647c6ae99SBarry Smith   PetscFunctionBegin;
1487aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
1488aa219208SBarry Smith   ierr = DMDASetDim(*da, 3);CHKERRQ(ierr);
1489aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, P);CHKERRQ(ierr);
1490aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, p);CHKERRQ(ierr);
14911321219cSEthan Coon   ierr = DMDASetBoundaryType(*da, bx, by, bz);CHKERRQ(ierr);
1492aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
1493aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
1494aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
1495aa219208SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, lz);CHKERRQ(ierr);
149647c6ae99SBarry Smith   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
14979a42bb27SBarry Smith   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
14989a42bb27SBarry Smith   ierr = DMSetUp(*da);CHKERRQ(ierr);
1499ca266f36SBarry Smith   ierr = DMViewFromOptions(*da,"-dm_view");CHKERRQ(ierr);
150047c6ae99SBarry Smith   PetscFunctionReturn(0);
150147c6ae99SBarry Smith }
1502