xref: /petsc/src/dm/impls/da/da3.c (revision e0877f539457ad1ce8178bc0750eac5ffa490a67)
17d0a6c19SBarry Smith 
247c6ae99SBarry Smith /*
347c6ae99SBarry Smith    Code for manipulating distributed regular 3d arrays in parallel.
447c6ae99SBarry Smith    File created by Peter Mell  7/14/95
547c6ae99SBarry Smith  */
647c6ae99SBarry Smith 
7af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h>     /*I   "petscdmda.h"    I*/
847c6ae99SBarry Smith 
99804daf3SBarry Smith #include <petscdraw.h>
1047c6ae99SBarry Smith #undef __FUNCT__
119a42bb27SBarry Smith #define __FUNCT__ "DMView_DA_3d"
12*e0877f53SBarry Smith static PetscErrorCode DMView_DA_3d(DM da,PetscViewer viewer)
1347c6ae99SBarry Smith {
1447c6ae99SBarry Smith   PetscErrorCode ierr;
1547c6ae99SBarry Smith   PetscMPIInt    rank;
169a42bb27SBarry Smith   PetscBool      iascii,isdraw,isbinary;
1747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
189a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
199a42bb27SBarry Smith   PetscBool ismatlab;
209a42bb27SBarry Smith #endif
2147c6ae99SBarry Smith 
2247c6ae99SBarry Smith   PetscFunctionBegin;
23ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr);
2447c6ae99SBarry Smith 
25251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
26251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
27251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
289a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
29251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr);
309a42bb27SBarry Smith #endif
3147c6ae99SBarry Smith   if (iascii) {
3247c6ae99SBarry Smith     PetscViewerFormat format;
3347c6ae99SBarry Smith 
341575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
3547c6ae99SBarry Smith     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
3647c6ae99SBarry Smith     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) {
37aa219208SBarry Smith       DMDALocalInfo info;
38aa219208SBarry Smith       ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
3947c6ae99SBarry 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);
4047c6ae99SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D, Y range of indices: %D %D, Z range of indices: %D %D\n",
4147c6ae99SBarry Smith                                                 info.xs,info.xs+info.xm,info.ys,info.ys+info.ym,info.zs,info.zs+info.zm);CHKERRQ(ierr);
4247c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX)
436636e97aSMatthew G Knepley       if (da->coordinates) {
4447c6ae99SBarry Smith         PetscInt        last;
4547c6ae99SBarry Smith         const PetscReal *coors;
466636e97aSMatthew G Knepley         ierr = VecGetArrayRead(da->coordinates,&coors);CHKERRQ(ierr);
476636e97aSMatthew G Knepley         ierr = VecGetLocalSize(da->coordinates,&last);CHKERRQ(ierr);
4847c6ae99SBarry Smith         last = last - 3;
4957622a8eSBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Lower left corner %g %g %g : Upper right %g %g %g\n",(double)coors[0],(double)coors[1],(double)coors[2],(double)coors[last],(double)coors[last+1],(double)coors[last+2]);CHKERRQ(ierr);
506636e97aSMatthew G Knepley         ierr = VecRestoreArrayRead(da->coordinates,&coors);CHKERRQ(ierr);
5147c6ae99SBarry Smith       }
5247c6ae99SBarry Smith #endif
5347c6ae99SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
541575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
55616157d6SJed Brown     } else {
56616157d6SJed Brown       ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr);
5747c6ae99SBarry Smith     }
5847c6ae99SBarry Smith   } else if (isdraw) {
5947c6ae99SBarry Smith     PetscDraw      draw;
6047c6ae99SBarry Smith     PetscReal      ymin = -1.0,ymax = (PetscReal)dd->N;
6147c6ae99SBarry Smith     PetscReal      xmin = -1.0,xmax = (PetscReal)((dd->M+2)*dd->P),x,y,ycoord,xcoord;
628ea3bf28SBarry Smith     PetscInt       k,plane,base;
638ea3bf28SBarry Smith     const PetscInt *idx;
6447c6ae99SBarry Smith     char           node[10];
6547c6ae99SBarry Smith     PetscBool      isnull;
6647c6ae99SBarry Smith 
6747c6ae99SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
6847c6ae99SBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
6947c6ae99SBarry Smith     ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
7047c6ae99SBarry Smith     ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
7147c6ae99SBarry Smith 
7247c6ae99SBarry Smith     /* first processor draw all node lines */
7347c6ae99SBarry Smith     if (!rank) {
7447c6ae99SBarry Smith       for (k=0; k<dd->P; k++) {
7547c6ae99SBarry Smith         ymin = 0.0; ymax = (PetscReal)(dd->N - 1);
7647c6ae99SBarry Smith         for (xmin=(PetscReal)(k*(dd->M+1)); xmin<(PetscReal)(dd->M+(k*(dd->M+1))); xmin++) {
7747c6ae99SBarry Smith           ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
7847c6ae99SBarry Smith         }
7947c6ae99SBarry Smith 
8047c6ae99SBarry Smith         xmin = (PetscReal)(k*(dd->M+1)); xmax = xmin + (PetscReal)(dd->M - 1);
8147c6ae99SBarry Smith         for (ymin=0; ymin<(PetscReal)dd->N; ymin++) {
8247c6ae99SBarry Smith           ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
8347c6ae99SBarry Smith         }
8447c6ae99SBarry Smith       }
8547c6ae99SBarry Smith     }
8647c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
8747c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
8847c6ae99SBarry Smith 
8947c6ae99SBarry Smith     for (k=0; k<dd->P; k++) {  /*Go through and draw for each plane*/
9047c6ae99SBarry Smith       if ((k >= dd->zs) && (k < dd->ze)) {
9147c6ae99SBarry Smith         /* draw my box */
9247c6ae99SBarry Smith         ymin = dd->ys;
9347c6ae99SBarry Smith         ymax = dd->ye - 1;
9447c6ae99SBarry Smith         xmin = dd->xs/dd->w    + (dd->M+1)*k;
9547c6ae99SBarry Smith         xmax =(dd->xe-1)/dd->w + (dd->M+1)*k;
9647c6ae99SBarry Smith 
9747c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
9847c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
9947c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
10047c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
10147c6ae99SBarry Smith 
10247c6ae99SBarry Smith         xmin = dd->xs/dd->w;
10347c6ae99SBarry Smith         xmax =(dd->xe-1)/dd->w;
10447c6ae99SBarry Smith 
10547c6ae99SBarry Smith         /* put in numbers*/
10647c6ae99SBarry Smith         base = (dd->base+(dd->xe-dd->xs)*(dd->ye-dd->ys)*(k-dd->zs))/dd->w;
10747c6ae99SBarry Smith 
10847c6ae99SBarry Smith         /* Identify which processor owns the box */
10947c6ae99SBarry Smith         sprintf(node,"%d",rank);
11047c6ae99SBarry Smith         ierr = PetscDrawString(draw,xmin+(dd->M+1)*k+.2,ymin+.3,PETSC_DRAW_RED,node);CHKERRQ(ierr);
11147c6ae99SBarry Smith 
11247c6ae99SBarry Smith         for (y=ymin; y<=ymax; y++) {
11347c6ae99SBarry Smith           for (x=xmin+(dd->M+1)*k; x<=xmax+(dd->M+1)*k; x++) {
11447c6ae99SBarry Smith             sprintf(node,"%d",(int)base++);
11547c6ae99SBarry Smith             ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
11647c6ae99SBarry Smith           }
11747c6ae99SBarry Smith         }
11847c6ae99SBarry Smith 
11947c6ae99SBarry Smith       }
12047c6ae99SBarry Smith     }
12147c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
12247c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
12347c6ae99SBarry Smith 
12447c6ae99SBarry Smith     for (k=0-dd->s; k<dd->P+dd->s; k++) {
12547c6ae99SBarry Smith       /* Go through and draw for each plane */
12647c6ae99SBarry Smith       if ((k >= dd->Zs) && (k < dd->Ze)) {
12747c6ae99SBarry Smith 
12847c6ae99SBarry Smith         /* overlay ghost numbers, useful for error checking */
129565245c5SBarry Smith         base = (dd->Xe-dd->Xs)*(dd->Ye-dd->Ys)*(k-dd->Zs)/dd->w;
13045b6f7e9SBarry Smith         ierr = ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr);
13147c6ae99SBarry Smith         plane=k;
132565245c5SBarry Smith         /* Keep z wrap around points on the drawing */
1338865f1eaSKarl Rupp         if (k<0) plane=dd->P+k;
1348865f1eaSKarl Rupp         if (k>=dd->P) plane=k-dd->P;
13547c6ae99SBarry Smith         ymin = dd->Ys; ymax = dd->Ye;
13647c6ae99SBarry Smith         xmin = (dd->M+1)*plane*dd->w;
13747c6ae99SBarry Smith         xmax = (dd->M+1)*plane*dd->w+dd->M*dd->w;
13847c6ae99SBarry Smith         for (y=ymin; y<ymax; y++) {
13947c6ae99SBarry Smith           for (x=xmin+dd->Xs; x<xmin+dd->Xe; x+=dd->w) {
1408ea3bf28SBarry Smith             sprintf(node,"%d",(int)(idx[base]));
14147c6ae99SBarry Smith             ycoord = y;
14247c6ae99SBarry Smith             /*Keep y wrap around points on drawing */
1438865f1eaSKarl Rupp             if (y<0) ycoord = dd->N+y;
14447c6ae99SBarry Smith 
1458865f1eaSKarl Rupp             if (y>=dd->N) ycoord = y-dd->N;
14647c6ae99SBarry Smith             xcoord = x;   /* Keep x wrap points on drawing */
14747c6ae99SBarry Smith 
1488865f1eaSKarl Rupp             if (x<xmin) xcoord = xmax - (xmin-x);
1498865f1eaSKarl Rupp             if (x>=xmax) xcoord = xmin + (x-xmax);
15047c6ae99SBarry Smith             ierr = PetscDrawString(draw,xcoord/dd->w,ycoord,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
151565245c5SBarry Smith             base++;
15247c6ae99SBarry Smith           }
15347c6ae99SBarry Smith         }
15445b6f7e9SBarry Smith         ierr = ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr);
15547c6ae99SBarry Smith       }
15647c6ae99SBarry Smith     }
15747c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
15847c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
1599a42bb27SBarry Smith   } else if (isbinary) {
1609a42bb27SBarry Smith     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
1619a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1629a42bb27SBarry Smith   } else if (ismatlab) {
1639a42bb27SBarry Smith     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
1649a42bb27SBarry Smith #endif
16511aeaf0aSBarry Smith   }
16647c6ae99SBarry Smith   PetscFunctionReturn(0);
16747c6ae99SBarry Smith }
16847c6ae99SBarry Smith 
16947c6ae99SBarry Smith #undef __FUNCT__
1709a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_3D"
1717087cfbeSBarry Smith PetscErrorCode  DMSetUp_DA_3D(DM da)
17247c6ae99SBarry Smith {
17347c6ae99SBarry Smith   DM_DA            *dd          = (DM_DA*)da->data;
17447c6ae99SBarry Smith   const PetscInt   M            = dd->M;
17547c6ae99SBarry Smith   const PetscInt   N            = dd->N;
17647c6ae99SBarry Smith   const PetscInt   P            = dd->P;
17747c6ae99SBarry Smith   PetscInt         m            = dd->m;
17847c6ae99SBarry Smith   PetscInt         n            = dd->n;
17947c6ae99SBarry Smith   PetscInt         p            = dd->p;
18047c6ae99SBarry Smith   const PetscInt   dof          = dd->w;
18147c6ae99SBarry Smith   const PetscInt   s            = dd->s;
182bff4a2f0SMatthew G. Knepley   DMBoundaryType   bx           = dd->bx;
183bff4a2f0SMatthew G. Knepley   DMBoundaryType   by           = dd->by;
184bff4a2f0SMatthew G. Knepley   DMBoundaryType   bz           = dd->bz;
18519fd82e9SBarry Smith   DMDAStencilType  stencil_type = dd->stencil_type;
18647c6ae99SBarry Smith   PetscInt         *lx          = dd->lx;
18747c6ae99SBarry Smith   PetscInt         *ly          = dd->ly;
18847c6ae99SBarry Smith   PetscInt         *lz          = dd->lz;
18947c6ae99SBarry Smith   MPI_Comm         comm;
19047c6ae99SBarry Smith   PetscMPIInt      rank,size;
191ce00eea3SSatish Balay   PetscInt         xs = 0,xe,ys = 0,ye,zs = 0,ze,x = 0,y = 0,z = 0;
192bd1fc5aeSBarry Smith   PetscInt         Xs,Xe,Ys,Ye,Zs,Ze,IXs,IXe,IYs,IYe,IZs,IZe,pm;
1938ea3bf28SBarry Smith   PetscInt         left,right,up,down,bottom,top,i,j,k,*idx,nn;
19447c6ae99SBarry Smith   PetscInt         n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14;
19547c6ae99SBarry Smith   PetscInt         n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26;
196db87c5ecSEthan Coon   PetscInt         *bases,*ldims,base,x_t,y_t,z_t,s_t,count,s_x,s_y,s_z;
19747c6ae99SBarry Smith   PetscInt         sn0  = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0;
19847c6ae99SBarry Smith   PetscInt         sn8  = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0;
19947c6ae99SBarry Smith   PetscInt         sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0;
20047c6ae99SBarry Smith   Vec              local,global;
201bd1fc5aeSBarry Smith   VecScatter       gtol;
20245b6f7e9SBarry Smith   IS               to,from;
2036f951b95Secoon   PetscBool        twod;
20447c6ae99SBarry Smith   PetscErrorCode   ierr;
20547c6ae99SBarry Smith 
2066f951b95Secoon 
20747c6ae99SBarry Smith   PetscFunctionBegin;
208bff4a2f0SMatthew G. Knepley   if (stencil_type == DMDA_STENCIL_BOX && (bx == DM_BOUNDARY_MIRROR || by == DM_BOUNDARY_MIRROR || bz == DM_BOUNDARY_MIRROR)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Mirror boundary and box stencil");
20947c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr);
2103855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES)
2113855c12bSBarry 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);
2123855c12bSBarry Smith #endif
2133855c12bSBarry Smith 
21447c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
21547c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
21647c6ae99SBarry Smith 
21747c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
21847c6ae99SBarry Smith     if (m < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
21947c6ae99SBarry Smith     else if (m > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
22047c6ae99SBarry Smith   }
22147c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
22247c6ae99SBarry Smith     if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
22347c6ae99SBarry Smith     else if (n > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
22447c6ae99SBarry Smith   }
22547c6ae99SBarry Smith   if (p != PETSC_DECIDE) {
22647c6ae99SBarry Smith     if (p < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %D",p);
22747c6ae99SBarry Smith     else if (p > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %D %d",p,size);
22847c6ae99SBarry Smith   }
2290716a85fSBarry 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);
23047c6ae99SBarry Smith 
23147c6ae99SBarry Smith   /* Partition the array among the processors */
23247c6ae99SBarry Smith   if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) {
23347c6ae99SBarry Smith     m = size/(n*p);
23447c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
23547c6ae99SBarry Smith     n = size/(m*p);
23647c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
23747c6ae99SBarry Smith     p = size/(m*n);
23847c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
23947c6ae99SBarry Smith     /* try for squarish distribution */
240369cc0aeSBarry Smith     m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N*p)));
24147c6ae99SBarry Smith     if (!m) m = 1;
24247c6ae99SBarry Smith     while (m > 0) {
24347c6ae99SBarry Smith       n = size/(m*p);
24447c6ae99SBarry Smith       if (m*n*p == size) break;
24547c6ae99SBarry Smith       m--;
24647c6ae99SBarry Smith     }
24747c6ae99SBarry Smith     if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %D",p);
24847c6ae99SBarry Smith     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
24947c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
25047c6ae99SBarry Smith     /* try for squarish distribution */
251369cc0aeSBarry Smith     m = (int)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n)));
25247c6ae99SBarry Smith     if (!m) m = 1;
25347c6ae99SBarry Smith     while (m > 0) {
25447c6ae99SBarry Smith       p = size/(m*n);
25547c6ae99SBarry Smith       if (m*n*p == size) break;
25647c6ae99SBarry Smith       m--;
25747c6ae99SBarry Smith     }
25847c6ae99SBarry Smith     if (!m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %D",n);
25947c6ae99SBarry Smith     if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;}
26047c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
26147c6ae99SBarry Smith     /* try for squarish distribution */
262369cc0aeSBarry Smith     n = (int)(0.5 + PetscSqrtReal(((PetscReal)N)*((PetscReal)size)/((PetscReal)P*m)));
26347c6ae99SBarry Smith     if (!n) n = 1;
26447c6ae99SBarry Smith     while (n > 0) {
26547c6ae99SBarry Smith       p = size/(m*n);
26647c6ae99SBarry Smith       if (m*n*p == size) break;
26747c6ae99SBarry Smith       n--;
26847c6ae99SBarry Smith     }
26947c6ae99SBarry Smith     if (!n) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %D",n);
27047c6ae99SBarry Smith     if (N > P && n < p) {PetscInt _n = n; n = p; p = _n;}
27147c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
27247c6ae99SBarry Smith     /* try for squarish distribution */
273369cc0aeSBarry Smith     n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),(PetscReal)(1./3.)));
27447c6ae99SBarry Smith     if (!n) n = 1;
27547c6ae99SBarry Smith     while (n > 0) {
27647c6ae99SBarry Smith       pm = size/n;
27747c6ae99SBarry Smith       if (n*pm == size) break;
27847c6ae99SBarry Smith       n--;
27947c6ae99SBarry Smith     }
28047c6ae99SBarry Smith     if (!n) n = 1;
281369cc0aeSBarry Smith     m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n)));
28247c6ae99SBarry Smith     if (!m) m = 1;
28347c6ae99SBarry Smith     while (m > 0) {
28447c6ae99SBarry Smith       p = size/(m*n);
28547c6ae99SBarry Smith       if (m*n*p == size) break;
28647c6ae99SBarry Smith       m--;
28747c6ae99SBarry Smith     }
28847c6ae99SBarry Smith     if (M > P && m < p) {PetscInt _m = m; m = p; p = _m;}
289ce94432eSBarry Smith   } else if (m*n*p != size) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
29047c6ae99SBarry Smith 
291ce94432eSBarry Smith   if (m*n*p != size) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_PLIB,"Could not find good partition");
292ce94432eSBarry Smith   if (M < m) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
293ce94432eSBarry Smith   if (N < n) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
294ce94432eSBarry Smith   if (P < p) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %D %D",P,p);
29547c6ae99SBarry Smith 
29647c6ae99SBarry Smith   /*
29747c6ae99SBarry Smith      Determine locally owned region
29847c6ae99SBarry Smith      [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes
29947c6ae99SBarry Smith   */
30047c6ae99SBarry Smith 
30147c6ae99SBarry Smith   if (!lx) {
302785e854fSJed Brown     ierr = PetscMalloc1(m, &dd->lx);CHKERRQ(ierr);
30347c6ae99SBarry Smith     lx   = dd->lx;
3048865f1eaSKarl Rupp     for (i=0; i<m; i++) lx[i] = M/m + ((M % m) > (i % m));
30547c6ae99SBarry Smith   }
30647c6ae99SBarry Smith   x  = lx[rank % m];
30747c6ae99SBarry Smith   xs = 0;
3088865f1eaSKarl Rupp   for (i=0; i<(rank%m); i++) xs += lx[i];
309bff4a2f0SMatthew G. Knepley   if ((x < s) && ((m > 1) || (bx == DM_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s);
31047c6ae99SBarry Smith 
31147c6ae99SBarry Smith   if (!ly) {
312785e854fSJed Brown     ierr = PetscMalloc1(n, &dd->ly);CHKERRQ(ierr);
31347c6ae99SBarry Smith     ly   = dd->ly;
3148865f1eaSKarl Rupp     for (i=0; i<n; i++) ly[i] = N/n + ((N % n) > (i % n));
31547c6ae99SBarry Smith   }
31647c6ae99SBarry Smith   y = ly[(rank % (m*n))/m];
317bff4a2f0SMatthew G. Knepley   if ((y < s) && ((n > 1) || (by == DM_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s);
31830729d88SBarry Smith 
31947c6ae99SBarry Smith   ys = 0;
3208865f1eaSKarl Rupp   for (i=0; i<(rank % (m*n))/m; i++) ys += ly[i];
32147c6ae99SBarry Smith 
32247c6ae99SBarry Smith   if (!lz) {
323785e854fSJed Brown     ierr = PetscMalloc1(p, &dd->lz);CHKERRQ(ierr);
32447c6ae99SBarry Smith     lz = dd->lz;
3258865f1eaSKarl Rupp     for (i=0; i<p; i++) lz[i] = P/p + ((P % p) > (i % p));
32647c6ae99SBarry Smith   }
32747c6ae99SBarry Smith   z = lz[rank/(m*n)];
328bcea557cSEthan Coon 
329fdc81ce8SEthan Coon   /* note this is different than x- and y-, as we will handle as an important special
330bff4a2f0SMatthew G. Knepley    case when p=P=1 and DM_BOUNDARY_PERIODIC and s > z.  This is to deal with 2D problems
331fdc81ce8SEthan Coon    in a 3D code.  Additional code for this case is noted with "2d case" comments */
3326f951b95Secoon   twod = PETSC_FALSE;
3338865f1eaSKarl Rupp   if (P == 1) twod = PETSC_TRUE;
334bff4a2f0SMatthew G. Knepley   else if ((z < s) && ((p > 1) || (bz == DM_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local z-width of domain z %D is smaller than stencil width s %D",z,s);
33547c6ae99SBarry Smith   zs = 0;
3368865f1eaSKarl Rupp   for (i=0; i<(rank/(m*n)); i++) zs += lz[i];
33747c6ae99SBarry Smith   ye = ys + y;
33847c6ae99SBarry Smith   xe = xs + x;
33947c6ae99SBarry Smith   ze = zs + z;
34047c6ae99SBarry Smith 
34188661749SPeter Brune   /* determine ghost region (Xs) and region scattered into (IXs)  */
342d9c9ebe5SPeter Brune   if (xs-s > 0) {
343d9c9ebe5SPeter Brune     Xs = xs - s; IXs = xs - s;
34488661749SPeter Brune   } else {
3458865f1eaSKarl Rupp     if (bx) Xs = xs - s;
3468865f1eaSKarl Rupp     else Xs = 0;
34788661749SPeter Brune     IXs = 0;
34888661749SPeter Brune   }
349d9c9ebe5SPeter Brune   if (xe+s <= M) {
350d9c9ebe5SPeter Brune     Xe = xe + s; IXe = xe + s;
35188661749SPeter Brune   } else {
35288661749SPeter Brune     if (bx) {
353d9c9ebe5SPeter Brune       Xs = xs - s; Xe = xe + s;
3548865f1eaSKarl Rupp     } else Xe = M;
35588661749SPeter Brune     IXe = M;
35688661749SPeter Brune   }
35747c6ae99SBarry Smith 
358bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) {
359d9c9ebe5SPeter Brune     IXs = xs - s;
360d9c9ebe5SPeter Brune     IXe = xe + s;
361d9c9ebe5SPeter Brune     Xs  = xs - s;
362d9c9ebe5SPeter Brune     Xe  = xe + s;
36388661749SPeter Brune   }
36488661749SPeter Brune 
365d9c9ebe5SPeter Brune   if (ys-s > 0) {
366d9c9ebe5SPeter Brune     Ys = ys - s; IYs = ys - s;
36788661749SPeter Brune   } else {
3688865f1eaSKarl Rupp     if (by) Ys = ys - s;
3698865f1eaSKarl Rupp     else Ys = 0;
37088661749SPeter Brune     IYs = 0;
37188661749SPeter Brune   }
372d9c9ebe5SPeter Brune   if (ye+s <= N) {
373d9c9ebe5SPeter Brune     Ye = ye + s; IYe = ye + s;
37488661749SPeter Brune   } else {
3758865f1eaSKarl Rupp     if (by) Ye = ye + s;
3768865f1eaSKarl Rupp     else Ye = N;
37788661749SPeter Brune     IYe = N;
37888661749SPeter Brune   }
37988661749SPeter Brune 
380bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) {
381d9c9ebe5SPeter Brune     IYs = ys - s;
382d9c9ebe5SPeter Brune     IYe = ye + s;
383d9c9ebe5SPeter Brune     Ys  = ys - s;
384d9c9ebe5SPeter Brune     Ye  = ye + s;
38588661749SPeter Brune   }
38688661749SPeter Brune 
387d9c9ebe5SPeter Brune   if (zs-s > 0) {
388d9c9ebe5SPeter Brune     Zs = zs - s; IZs = zs - s;
38988661749SPeter Brune   } else {
3908865f1eaSKarl Rupp     if (bz) Zs = zs - s;
3918865f1eaSKarl Rupp     else Zs = 0;
39288661749SPeter Brune     IZs = 0;
39388661749SPeter Brune   }
394d9c9ebe5SPeter Brune   if (ze+s <= P) {
395d9c9ebe5SPeter Brune     Ze = ze + s; IZe = ze + s;
39688661749SPeter Brune   } else {
3978865f1eaSKarl Rupp     if (bz) Ze = ze + s;
3988865f1eaSKarl Rupp     else Ze = P;
39988661749SPeter Brune     IZe = P;
40088661749SPeter Brune   }
40188661749SPeter Brune 
402bff4a2f0SMatthew G. Knepley   if (bz == DM_BOUNDARY_PERIODIC || bz == DM_BOUNDARY_MIRROR) {
403d9c9ebe5SPeter Brune     IZs = zs - s;
404d9c9ebe5SPeter Brune     IZe = ze + s;
405d9c9ebe5SPeter Brune     Zs  = zs - s;
406d9c9ebe5SPeter Brune     Ze  = ze + s;
40788661749SPeter Brune   }
40847c6ae99SBarry Smith 
40947c6ae99SBarry Smith   /* Resize all X parameters to reflect w */
410d9c9ebe5SPeter Brune   s_x = s;
411d9c9ebe5SPeter Brune   s_y = s;
412d9c9ebe5SPeter Brune   s_z = s;
41347c6ae99SBarry Smith 
41447c6ae99SBarry Smith   /* determine starting point of each processor */
41547c6ae99SBarry Smith   nn       = x*y*z;
416dcca6d9dSJed Brown   ierr     = PetscMalloc2(size+1,&bases,size,&ldims);CHKERRQ(ierr);
41747c6ae99SBarry Smith   ierr     = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
41847c6ae99SBarry Smith   bases[0] = 0;
4198865f1eaSKarl Rupp   for (i=1; i<=size; i++) bases[i] = ldims[i-1];
4208865f1eaSKarl Rupp   for (i=1; i<=size; i++) bases[i] += bases[i-1];
421ce00eea3SSatish Balay   base = bases[rank]*dof;
42247c6ae99SBarry Smith 
42347c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
424ce00eea3SSatish Balay   dd->Nlocal = x*y*z*dof;
425b1fb7eb7SBarry Smith   ierr       = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global);CHKERRQ(ierr);
426ce00eea3SSatish Balay   dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs)*dof;
427b1fb7eb7SBarry Smith   ierr       = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local);CHKERRQ(ierr);
42847c6ae99SBarry Smith 
42947c6ae99SBarry Smith   /* generate appropriate vector scatters */
43047c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
43102fe608eSBarry Smith   ierr   = PetscMalloc1((IXe-IXs)*(IYe-IYs)*(IZe-IZs),&idx);CHKERRQ(ierr);
432ce00eea3SSatish Balay   left   = xs - Xs; right = left + x;
43347c6ae99SBarry Smith   bottom = ys - Ys; top = bottom + y;
43447c6ae99SBarry Smith   down   = zs - Zs; up  = down + z;
43547c6ae99SBarry Smith   count  = 0;
43647c6ae99SBarry Smith   for (i=down; i<up; i++) {
43747c6ae99SBarry Smith     for (j=bottom; j<top; j++) {
438ce00eea3SSatish Balay       for (k=left; k<right; k++) {
439ce00eea3SSatish Balay         idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
44047c6ae99SBarry Smith       }
44147c6ae99SBarry Smith     }
44247c6ae99SBarry Smith   }
44347c6ae99SBarry Smith 
444ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
445ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
446d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_BOX) {
447ce00eea3SSatish Balay     left   = IXs - Xs; right = left + (IXe-IXs);
448ce00eea3SSatish Balay     bottom = IYs - Ys; top = bottom + (IYe-IYs);
449ce00eea3SSatish Balay     down   = IZs - Zs; up  = down + (IZe-IZs);
450ce00eea3SSatish Balay     count  = 0;
451ce00eea3SSatish Balay     for (i=down; i<up; i++) {
452ce00eea3SSatish Balay       for (j=bottom; j<top; j++) {
453ce00eea3SSatish Balay         for (k=left; k<right; k++) {
454ce00eea3SSatish Balay           idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
455ce00eea3SSatish Balay         }
456ce00eea3SSatish Balay       }
457ce00eea3SSatish Balay     }
458ce00eea3SSatish Balay     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
45947c6ae99SBarry Smith   } else {
46047c6ae99SBarry Smith     /* This is way ugly! We need to list the funny cross type region */
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;
465ce00eea3SSatish Balay     /* the bottom chunck */
466ce00eea3SSatish Balay     for (i=(IZs-Zs); i<down; i++) {
46747c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
468ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
46947c6ae99SBarry Smith       }
47047c6ae99SBarry Smith     }
47147c6ae99SBarry Smith     /* the middle piece */
47247c6ae99SBarry Smith     for (i=down; i<up; i++) {
47347c6ae99SBarry Smith       /* front */
474ce00eea3SSatish Balay       for (j=(IYs-Ys); j<bottom; j++) {
475ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
47647c6ae99SBarry Smith       }
47747c6ae99SBarry Smith       /* middle */
47847c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
479ce00eea3SSatish Balay         for (k=IXs-Xs; k<IXe-Xs; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
48047c6ae99SBarry Smith       }
48147c6ae99SBarry Smith       /* back */
482ce00eea3SSatish Balay       for (j=top; j<top+IYe-ye; j++) {
483ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
48447c6ae99SBarry Smith       }
48547c6ae99SBarry Smith     }
48647c6ae99SBarry Smith     /* the top piece */
487ce00eea3SSatish Balay     for (i=up; i<up+IZe-ze; i++) {
48847c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
489ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
49047c6ae99SBarry Smith       }
49147c6ae99SBarry Smith     }
49247c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
49347c6ae99SBarry Smith   }
49447c6ae99SBarry Smith 
49547c6ae99SBarry Smith   /* determine who lies on each side of use stored in    n24 n25 n26
49647c6ae99SBarry Smith                                                          n21 n22 n23
49747c6ae99SBarry Smith                                                          n18 n19 n20
49847c6ae99SBarry Smith 
49947c6ae99SBarry Smith                                                          n15 n16 n17
50047c6ae99SBarry Smith                                                          n12     n14
50147c6ae99SBarry Smith                                                          n9  n10 n11
50247c6ae99SBarry Smith 
50347c6ae99SBarry Smith                                                          n6  n7  n8
50447c6ae99SBarry Smith                                                          n3  n4  n5
50547c6ae99SBarry Smith                                                          n0  n1  n2
50647c6ae99SBarry Smith   */
50747c6ae99SBarry Smith 
50847c6ae99SBarry Smith   /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
50947c6ae99SBarry Smith   /* Assume Nodes are Internal to the Cube */
51047c6ae99SBarry Smith   n0 = rank - m*n - m - 1;
51147c6ae99SBarry Smith   n1 = rank - m*n - m;
51247c6ae99SBarry Smith   n2 = rank - m*n - m + 1;
51347c6ae99SBarry Smith   n3 = rank - m*n -1;
51447c6ae99SBarry Smith   n4 = rank - m*n;
51547c6ae99SBarry Smith   n5 = rank - m*n + 1;
51647c6ae99SBarry Smith   n6 = rank - m*n + m - 1;
51747c6ae99SBarry Smith   n7 = rank - m*n + m;
51847c6ae99SBarry Smith   n8 = rank - m*n + m + 1;
51947c6ae99SBarry Smith 
52047c6ae99SBarry Smith   n9  = rank - m - 1;
52147c6ae99SBarry Smith   n10 = rank - m;
52247c6ae99SBarry Smith   n11 = rank - m + 1;
52347c6ae99SBarry Smith   n12 = rank - 1;
52447c6ae99SBarry Smith   n14 = rank + 1;
52547c6ae99SBarry Smith   n15 = rank + m - 1;
52647c6ae99SBarry Smith   n16 = rank + m;
52747c6ae99SBarry Smith   n17 = rank + m + 1;
52847c6ae99SBarry Smith 
52947c6ae99SBarry Smith   n18 = rank + m*n - m - 1;
53047c6ae99SBarry Smith   n19 = rank + m*n - m;
53147c6ae99SBarry Smith   n20 = rank + m*n - m + 1;
53247c6ae99SBarry Smith   n21 = rank + m*n - 1;
53347c6ae99SBarry Smith   n22 = rank + m*n;
53447c6ae99SBarry Smith   n23 = rank + m*n + 1;
53547c6ae99SBarry Smith   n24 = rank + m*n + m - 1;
53647c6ae99SBarry Smith   n25 = rank + m*n + m;
53747c6ae99SBarry Smith   n26 = rank + m*n + m + 1;
53847c6ae99SBarry Smith 
53947c6ae99SBarry Smith   /* Assume Pieces are on Faces of Cube */
54047c6ae99SBarry Smith 
54147c6ae99SBarry Smith   if (xs == 0) { /* First assume not corner or edge */
54247c6ae99SBarry Smith     n0  = rank       -1 - (m*n);
54347c6ae99SBarry Smith     n3  = rank + m   -1 - (m*n);
54447c6ae99SBarry Smith     n6  = rank + 2*m -1 - (m*n);
54547c6ae99SBarry Smith     n9  = rank       -1;
54647c6ae99SBarry Smith     n12 = rank + m   -1;
54747c6ae99SBarry Smith     n15 = rank + 2*m -1;
54847c6ae99SBarry Smith     n18 = rank       -1 + (m*n);
54947c6ae99SBarry Smith     n21 = rank + m   -1 + (m*n);
55047c6ae99SBarry Smith     n24 = rank + 2*m -1 + (m*n);
55147c6ae99SBarry Smith   }
55247c6ae99SBarry Smith 
553ce00eea3SSatish Balay   if (xe == M) { /* First assume not corner or edge */
55447c6ae99SBarry Smith     n2  = rank -2*m +1 - (m*n);
55547c6ae99SBarry Smith     n5  = rank - m  +1 - (m*n);
55647c6ae99SBarry Smith     n8  = rank      +1 - (m*n);
55747c6ae99SBarry Smith     n11 = rank -2*m +1;
55847c6ae99SBarry Smith     n14 = rank - m  +1;
55947c6ae99SBarry Smith     n17 = rank      +1;
56047c6ae99SBarry Smith     n20 = rank -2*m +1 + (m*n);
56147c6ae99SBarry Smith     n23 = rank - m  +1 + (m*n);
56247c6ae99SBarry Smith     n26 = rank      +1 + (m*n);
56347c6ae99SBarry Smith   }
56447c6ae99SBarry Smith 
56547c6ae99SBarry Smith   if (ys==0) { /* First assume not corner or edge */
56647c6ae99SBarry Smith     n0  = rank + m * (n-1) -1 - (m*n);
56747c6ae99SBarry Smith     n1  = rank + m * (n-1)    - (m*n);
56847c6ae99SBarry Smith     n2  = rank + m * (n-1) +1 - (m*n);
56947c6ae99SBarry Smith     n9  = rank + m * (n-1) -1;
57047c6ae99SBarry Smith     n10 = rank + m * (n-1);
57147c6ae99SBarry Smith     n11 = rank + m * (n-1) +1;
57247c6ae99SBarry Smith     n18 = rank + m * (n-1) -1 + (m*n);
57347c6ae99SBarry Smith     n19 = rank + m * (n-1)    + (m*n);
57447c6ae99SBarry Smith     n20 = rank + m * (n-1) +1 + (m*n);
57547c6ae99SBarry Smith   }
57647c6ae99SBarry Smith 
57747c6ae99SBarry Smith   if (ye == N) { /* First assume not corner or edge */
57847c6ae99SBarry Smith     n6  = rank - m * (n-1) -1 - (m*n);
57947c6ae99SBarry Smith     n7  = rank - m * (n-1)    - (m*n);
58047c6ae99SBarry Smith     n8  = rank - m * (n-1) +1 - (m*n);
58147c6ae99SBarry Smith     n15 = rank - m * (n-1) -1;
58247c6ae99SBarry Smith     n16 = rank - m * (n-1);
58347c6ae99SBarry Smith     n17 = rank - m * (n-1) +1;
58447c6ae99SBarry Smith     n24 = rank - m * (n-1) -1 + (m*n);
58547c6ae99SBarry Smith     n25 = rank - m * (n-1)    + (m*n);
58647c6ae99SBarry Smith     n26 = rank - m * (n-1) +1 + (m*n);
58747c6ae99SBarry Smith   }
58847c6ae99SBarry Smith 
58947c6ae99SBarry Smith   if (zs == 0) { /* First assume not corner or edge */
59047c6ae99SBarry Smith     n0 = size - (m*n) + rank - m - 1;
59147c6ae99SBarry Smith     n1 = size - (m*n) + rank - m;
59247c6ae99SBarry Smith     n2 = size - (m*n) + rank - m + 1;
59347c6ae99SBarry Smith     n3 = size - (m*n) + rank - 1;
59447c6ae99SBarry Smith     n4 = size - (m*n) + rank;
59547c6ae99SBarry Smith     n5 = size - (m*n) + rank + 1;
59647c6ae99SBarry Smith     n6 = size - (m*n) + rank + m - 1;
59747c6ae99SBarry Smith     n7 = size - (m*n) + rank + m;
59847c6ae99SBarry Smith     n8 = size - (m*n) + rank + m + 1;
59947c6ae99SBarry Smith   }
60047c6ae99SBarry Smith 
60147c6ae99SBarry Smith   if (ze == P) { /* First assume not corner or edge */
60247c6ae99SBarry Smith     n18 = (m*n) - (size-rank) - m - 1;
60347c6ae99SBarry Smith     n19 = (m*n) - (size-rank) - m;
60447c6ae99SBarry Smith     n20 = (m*n) - (size-rank) - m + 1;
60547c6ae99SBarry Smith     n21 = (m*n) - (size-rank) - 1;
60647c6ae99SBarry Smith     n22 = (m*n) - (size-rank);
60747c6ae99SBarry Smith     n23 = (m*n) - (size-rank) + 1;
60847c6ae99SBarry Smith     n24 = (m*n) - (size-rank) + m - 1;
60947c6ae99SBarry Smith     n25 = (m*n) - (size-rank) + m;
61047c6ae99SBarry Smith     n26 = (m*n) - (size-rank) + m + 1;
61147c6ae99SBarry Smith   }
61247c6ae99SBarry Smith 
61347c6ae99SBarry Smith   if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */
61447c6ae99SBarry Smith     n0 = size - m*n + rank + m-1 - m;
61547c6ae99SBarry Smith     n3 = size - m*n + rank + m-1;
61647c6ae99SBarry Smith     n6 = size - m*n + rank + m-1 + m;
61747c6ae99SBarry Smith   }
61847c6ae99SBarry Smith 
61947c6ae99SBarry Smith   if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */
62047c6ae99SBarry Smith     n18 = m*n - (size - rank) + m-1 - m;
62147c6ae99SBarry Smith     n21 = m*n - (size - rank) + m-1;
62247c6ae99SBarry Smith     n24 = m*n - (size - rank) + m-1 + m;
62347c6ae99SBarry Smith   }
62447c6ae99SBarry Smith 
62547c6ae99SBarry Smith   if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */
62647c6ae99SBarry Smith     n0  = rank + m*n -1 - m*n;
62747c6ae99SBarry Smith     n9  = rank + m*n -1;
62847c6ae99SBarry Smith     n18 = rank + m*n -1 + m*n;
62947c6ae99SBarry Smith   }
63047c6ae99SBarry Smith 
63147c6ae99SBarry Smith   if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */
63247c6ae99SBarry Smith     n6  = rank - m*(n-1) + m-1 - m*n;
63347c6ae99SBarry Smith     n15 = rank - m*(n-1) + m-1;
63447c6ae99SBarry Smith     n24 = rank - m*(n-1) + m-1 + m*n;
63547c6ae99SBarry Smith   }
63647c6ae99SBarry Smith 
637ce00eea3SSatish Balay   if ((xe==M) && (zs==0)) { /* Assume an edge, not corner */
63847c6ae99SBarry Smith     n2 = size - (m*n-rank) - (m-1) - m;
63947c6ae99SBarry Smith     n5 = size - (m*n-rank) - (m-1);
64047c6ae99SBarry Smith     n8 = size - (m*n-rank) - (m-1) + m;
64147c6ae99SBarry Smith   }
64247c6ae99SBarry Smith 
643ce00eea3SSatish Balay   if ((xe==M) && (ze==P)) { /* Assume an edge, not corner */
64447c6ae99SBarry Smith     n20 = m*n - (size - rank) - (m-1) - m;
64547c6ae99SBarry Smith     n23 = m*n - (size - rank) - (m-1);
64647c6ae99SBarry Smith     n26 = m*n - (size - rank) - (m-1) + m;
64747c6ae99SBarry Smith   }
64847c6ae99SBarry Smith 
649ce00eea3SSatish Balay   if ((xe==M) && (ys==0)) { /* Assume an edge, not corner */
65047c6ae99SBarry Smith     n2  = rank + m*(n-1) - (m-1) - m*n;
65147c6ae99SBarry Smith     n11 = rank + m*(n-1) - (m-1);
65247c6ae99SBarry Smith     n20 = rank + m*(n-1) - (m-1) + m*n;
65347c6ae99SBarry Smith   }
65447c6ae99SBarry Smith 
655ce00eea3SSatish Balay   if ((xe==M) && (ye==N)) { /* Assume an edge, not corner */
65647c6ae99SBarry Smith     n8  = rank - m*n +1 - m*n;
65747c6ae99SBarry Smith     n17 = rank - m*n +1;
65847c6ae99SBarry Smith     n26 = rank - m*n +1 + m*n;
65947c6ae99SBarry Smith   }
66047c6ae99SBarry Smith 
66147c6ae99SBarry Smith   if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */
66247c6ae99SBarry Smith     n0 = size - m + rank -1;
66347c6ae99SBarry Smith     n1 = size - m + rank;
66447c6ae99SBarry Smith     n2 = size - m + rank +1;
66547c6ae99SBarry Smith   }
66647c6ae99SBarry Smith 
66747c6ae99SBarry Smith   if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */
66847c6ae99SBarry Smith     n18 = m*n - (size - rank) + m*(n-1) -1;
66947c6ae99SBarry Smith     n19 = m*n - (size - rank) + m*(n-1);
67047c6ae99SBarry Smith     n20 = m*n - (size - rank) + m*(n-1) +1;
67147c6ae99SBarry Smith   }
67247c6ae99SBarry Smith 
67347c6ae99SBarry Smith   if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */
67447c6ae99SBarry Smith     n6 = size - (m*n-rank) - m * (n-1) -1;
67547c6ae99SBarry Smith     n7 = size - (m*n-rank) - m * (n-1);
67647c6ae99SBarry Smith     n8 = size - (m*n-rank) - m * (n-1) +1;
67747c6ae99SBarry Smith   }
67847c6ae99SBarry Smith 
67947c6ae99SBarry Smith   if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */
68047c6ae99SBarry Smith     n24 = rank - (size-m) -1;
68147c6ae99SBarry Smith     n25 = rank - (size-m);
68247c6ae99SBarry Smith     n26 = rank - (size-m) +1;
68347c6ae99SBarry Smith   }
68447c6ae99SBarry Smith 
68547c6ae99SBarry Smith   /* Check for Corners */
6868865f1eaSKarl Rupp   if ((xs==0) && (ys==0) && (zs==0)) n0  = size -1;
6878865f1eaSKarl Rupp   if ((xs==0) && (ys==0) && (ze==P)) n18 = m*n-1;
6888865f1eaSKarl Rupp   if ((xs==0) && (ye==N) && (zs==0)) n6  = (size-1)-m*(n-1);
6898865f1eaSKarl Rupp   if ((xs==0) && (ye==N) && (ze==P)) n24 = m-1;
6908865f1eaSKarl Rupp   if ((xe==M) && (ys==0) && (zs==0)) n2  = size-m;
6918865f1eaSKarl Rupp   if ((xe==M) && (ys==0) && (ze==P)) n20 = m*n-m;
6928865f1eaSKarl Rupp   if ((xe==M) && (ye==N) && (zs==0)) n8  = size-m*n;
6938865f1eaSKarl Rupp   if ((xe==M) && (ye==N) && (ze==P)) n26 = 0;
69447c6ae99SBarry Smith 
69547c6ae99SBarry Smith   /* Check for when not X,Y, and Z Periodic */
69647c6ae99SBarry Smith 
69747c6ae99SBarry Smith   /* If not X periodic */
698bff4a2f0SMatthew G. Knepley   if (bx != DM_BOUNDARY_PERIODIC) {
6998865f1eaSKarl Rupp     if (xs==0) n0 = n3 = n6 = n9  = n12 = n15 = n18 = n21 = n24 = -2;
7008865f1eaSKarl Rupp     if (xe==M) n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2;
70147c6ae99SBarry Smith   }
70247c6ae99SBarry Smith 
70347c6ae99SBarry Smith   /* If not Y periodic */
704bff4a2f0SMatthew G. Knepley   if (by != DM_BOUNDARY_PERIODIC) {
7058865f1eaSKarl Rupp     if (ys==0) n0 = n1 = n2 = n9  = n10 = n11 = n18 = n19 = n20 = -2;
7068865f1eaSKarl Rupp     if (ye==N) n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2;
70747c6ae99SBarry Smith   }
70847c6ae99SBarry Smith 
70947c6ae99SBarry Smith   /* If not Z periodic */
710bff4a2f0SMatthew G. Knepley   if (bz != DM_BOUNDARY_PERIODIC) {
7118865f1eaSKarl Rupp     if (zs==0) n0  = n1  = n2  = n3  = n4  = n5  = n6  = n7  = n8  = -2;
7128865f1eaSKarl Rupp     if (ze==P) n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;
71347c6ae99SBarry Smith   }
71447c6ae99SBarry Smith 
715785e854fSJed Brown   ierr = PetscMalloc1(27,&dd->neighbors);CHKERRQ(ierr);
7168865f1eaSKarl Rupp 
71747c6ae99SBarry Smith   dd->neighbors[0]  = n0;
71847c6ae99SBarry Smith   dd->neighbors[1]  = n1;
71947c6ae99SBarry Smith   dd->neighbors[2]  = n2;
72047c6ae99SBarry Smith   dd->neighbors[3]  = n3;
72147c6ae99SBarry Smith   dd->neighbors[4]  = n4;
72247c6ae99SBarry Smith   dd->neighbors[5]  = n5;
72347c6ae99SBarry Smith   dd->neighbors[6]  = n6;
72447c6ae99SBarry Smith   dd->neighbors[7]  = n7;
72547c6ae99SBarry Smith   dd->neighbors[8]  = n8;
72647c6ae99SBarry Smith   dd->neighbors[9]  = n9;
72747c6ae99SBarry Smith   dd->neighbors[10] = n10;
72847c6ae99SBarry Smith   dd->neighbors[11] = n11;
72947c6ae99SBarry Smith   dd->neighbors[12] = n12;
73047c6ae99SBarry Smith   dd->neighbors[13] = rank;
73147c6ae99SBarry Smith   dd->neighbors[14] = n14;
73247c6ae99SBarry Smith   dd->neighbors[15] = n15;
73347c6ae99SBarry Smith   dd->neighbors[16] = n16;
73447c6ae99SBarry Smith   dd->neighbors[17] = n17;
73547c6ae99SBarry Smith   dd->neighbors[18] = n18;
73647c6ae99SBarry Smith   dd->neighbors[19] = n19;
73747c6ae99SBarry Smith   dd->neighbors[20] = n20;
73847c6ae99SBarry Smith   dd->neighbors[21] = n21;
73947c6ae99SBarry Smith   dd->neighbors[22] = n22;
74047c6ae99SBarry Smith   dd->neighbors[23] = n23;
74147c6ae99SBarry Smith   dd->neighbors[24] = n24;
74247c6ae99SBarry Smith   dd->neighbors[25] = n25;
74347c6ae99SBarry Smith   dd->neighbors[26] = n26;
74447c6ae99SBarry Smith 
74547c6ae99SBarry Smith   /* If star stencil then delete the corner neighbors */
746d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
74747c6ae99SBarry Smith     /* save information about corner neighbors */
74847c6ae99SBarry Smith     sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7;
74947c6ae99SBarry Smith     sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18;
75047c6ae99SBarry Smith     sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25;
75147c6ae99SBarry Smith     sn26 = n26;
7528865f1eaSKarl Rupp     n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 = n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1;
75347c6ae99SBarry Smith   }
75447c6ae99SBarry Smith 
755785e854fSJed Brown   ierr = PetscMalloc1((Xe-Xs)*(Ye-Ys)*(Ze-Zs),&idx);CHKERRQ(ierr);
75647c6ae99SBarry Smith 
75747c6ae99SBarry Smith   nn = 0;
75847c6ae99SBarry Smith   /* Bottom Level */
75947c6ae99SBarry Smith   for (k=0; k<s_z; k++) {
76047c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
76147c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
762ce00eea3SSatish Balay         x_t = lx[n0 % m];
76347c6ae99SBarry Smith         y_t = ly[(n0 % (m*n))/m];
76447c6ae99SBarry Smith         z_t = lz[n0 / (m*n)];
76547c6ae99SBarry 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;
7668865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n0] + x_t*y_t*z_t - (s_y-i)*x_t - s_x; /* 2D case */
7678865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
76847c6ae99SBarry Smith       }
76947c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
77047c6ae99SBarry Smith         x_t = x;
77147c6ae99SBarry Smith         y_t = ly[(n1 % (m*n))/m];
77247c6ae99SBarry Smith         z_t = lz[n1 / (m*n)];
77347c6ae99SBarry 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;
7748865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */
7758865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
77647c6ae99SBarry Smith       }
77747c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
778ce00eea3SSatish Balay         x_t = lx[n2 % m];
77947c6ae99SBarry Smith         y_t = ly[(n2 % (m*n))/m];
78047c6ae99SBarry Smith         z_t = lz[n2 / (m*n)];
78147c6ae99SBarry 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;
7828865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */
7838865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
78447c6ae99SBarry Smith       }
78547c6ae99SBarry Smith     }
78647c6ae99SBarry Smith 
78747c6ae99SBarry Smith     for (i=0; i<y; i++) {
78847c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
789ce00eea3SSatish Balay         x_t = lx[n3 % m];
79047c6ae99SBarry Smith         y_t = y;
79147c6ae99SBarry Smith         z_t = lz[n3 / (m*n)];
79247c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
7938865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - x_t*y_t; /* 2D case */
7948865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
79547c6ae99SBarry Smith       }
79647c6ae99SBarry Smith 
79747c6ae99SBarry Smith       if (n4 >= 0) { /* middle */
79847c6ae99SBarry Smith         x_t = x;
79947c6ae99SBarry Smith         y_t = y;
80047c6ae99SBarry Smith         z_t = lz[n4 / (m*n)];
80147c6ae99SBarry Smith         s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8028865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n4] + i*x_t + x_t*y_t*z_t - x_t*y_t; /* 2D case */
8038865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
804bff4a2f0SMatthew G. Knepley       } else if (bz == DM_BOUNDARY_MIRROR) {
8058865f1eaSKarl Rupp         for (j=0; j<x; j++) idx[nn++] = 0;
80647c6ae99SBarry Smith       }
80747c6ae99SBarry Smith 
80847c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
809ce00eea3SSatish Balay         x_t = lx[n5 % m];
81047c6ae99SBarry Smith         y_t = y;
81147c6ae99SBarry Smith         z_t = lz[n5 / (m*n)];
81247c6ae99SBarry Smith         s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8138865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n5] + i*x_t + x_t*y_t*z_t - x_t*y_t; /* 2D case */
8148865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
81547c6ae99SBarry Smith       }
81647c6ae99SBarry Smith     }
81747c6ae99SBarry Smith 
81847c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
81947c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
820ce00eea3SSatish Balay         x_t = lx[n6 % m];
82147c6ae99SBarry Smith         y_t = ly[(n6 % (m*n))/m];
82247c6ae99SBarry Smith         z_t = lz[n6 / (m*n)];
82347c6ae99SBarry Smith         s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8248865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - x_t*y_t; /* 2D case */
8258865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
82647c6ae99SBarry Smith       }
82747c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
82847c6ae99SBarry Smith         x_t = x;
82947c6ae99SBarry Smith         y_t = ly[(n7 % (m*n))/m];
83047c6ae99SBarry Smith         z_t = lz[n7 / (m*n)];
83147c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8328865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - x_t*y_t; /* 2D case */
8338865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
83447c6ae99SBarry Smith       }
83547c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
836ce00eea3SSatish Balay         x_t = lx[n8 % m];
83747c6ae99SBarry Smith         y_t = ly[(n8 % (m*n))/m];
83847c6ae99SBarry Smith         z_t = lz[n8 / (m*n)];
83947c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8408865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - x_t*y_t; /* 2D case */
8418865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
84247c6ae99SBarry Smith       }
84347c6ae99SBarry Smith     }
84447c6ae99SBarry Smith   }
84547c6ae99SBarry Smith 
84647c6ae99SBarry Smith   /* Middle Level */
84747c6ae99SBarry Smith   for (k=0; k<z; k++) {
84847c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
84947c6ae99SBarry Smith       if (n9 >= 0) { /* left below */
850ce00eea3SSatish Balay         x_t = lx[n9 % m];
85147c6ae99SBarry Smith         y_t = ly[(n9 % (m*n))/m];
85247c6ae99SBarry Smith         /* z_t = z; */
85347c6ae99SBarry Smith         s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
8548865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
85547c6ae99SBarry Smith       }
85647c6ae99SBarry Smith       if (n10 >= 0) { /* directly below */
85747c6ae99SBarry Smith         x_t = x;
85847c6ae99SBarry Smith         y_t = ly[(n10 % (m*n))/m];
85947c6ae99SBarry Smith         /* z_t = z; */
86047c6ae99SBarry Smith         s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
8618865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
862bff4a2f0SMatthew G. Knepley       }  else if (by == DM_BOUNDARY_MIRROR) {
8638865f1eaSKarl Rupp         for (j=0; j<x; j++) idx[nn++] = 0;
86447c6ae99SBarry Smith       }
86547c6ae99SBarry Smith       if (n11 >= 0) { /* right below */
866ce00eea3SSatish Balay         x_t = lx[n11 % m];
86747c6ae99SBarry Smith         y_t = ly[(n11 % (m*n))/m];
86847c6ae99SBarry Smith         /* z_t = z; */
86947c6ae99SBarry Smith         s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
8708865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
87147c6ae99SBarry Smith       }
87247c6ae99SBarry Smith     }
87347c6ae99SBarry Smith 
87447c6ae99SBarry Smith     for (i=0; i<y; i++) {
87547c6ae99SBarry Smith       if (n12 >= 0) { /* directly left */
876ce00eea3SSatish Balay         x_t = lx[n12 % m];
87747c6ae99SBarry Smith         y_t = y;
87847c6ae99SBarry Smith         /* z_t = z; */
87947c6ae99SBarry Smith         s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
8808865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
881bff4a2f0SMatthew G. Knepley       }  else if (bx == DM_BOUNDARY_MIRROR) {
8828865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = 0;
88347c6ae99SBarry Smith       }
88447c6ae99SBarry Smith 
88547c6ae99SBarry Smith       /* Interior */
88647c6ae99SBarry Smith       s_t = bases[rank] + i*x + k*x*y;
8878865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = s_t++;
88847c6ae99SBarry Smith 
88947c6ae99SBarry Smith       if (n14 >= 0) { /* directly right */
890ce00eea3SSatish Balay         x_t = lx[n14 % m];
89147c6ae99SBarry Smith         y_t = y;
89247c6ae99SBarry Smith         /* z_t = z; */
89347c6ae99SBarry Smith         s_t = bases[n14] + i*x_t + k*x_t*y_t;
8948865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
895bff4a2f0SMatthew G. Knepley       } else if (bx == DM_BOUNDARY_MIRROR) {
8968865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = 0;
89747c6ae99SBarry Smith       }
89847c6ae99SBarry Smith     }
89947c6ae99SBarry Smith 
90047c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
90147c6ae99SBarry Smith       if (n15 >= 0) { /* left above */
902ce00eea3SSatish Balay         x_t = lx[n15 % m];
90347c6ae99SBarry Smith         y_t = ly[(n15 % (m*n))/m];
90447c6ae99SBarry Smith         /* z_t = z; */
90547c6ae99SBarry Smith         s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
9068865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
90747c6ae99SBarry Smith       }
90847c6ae99SBarry Smith       if (n16 >= 0) { /* directly above */
90947c6ae99SBarry Smith         x_t = x;
91047c6ae99SBarry Smith         y_t = ly[(n16 % (m*n))/m];
91147c6ae99SBarry Smith         /* z_t = z; */
91247c6ae99SBarry Smith         s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
9138865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
914bff4a2f0SMatthew G. Knepley       } else if (by == DM_BOUNDARY_MIRROR) {
9158865f1eaSKarl Rupp         for (j=0; j<x; j++) idx[nn++] = 0;
91647c6ae99SBarry Smith       }
91747c6ae99SBarry Smith       if (n17 >= 0) { /* right above */
918ce00eea3SSatish Balay         x_t = lx[n17 % m];
91947c6ae99SBarry Smith         y_t = ly[(n17 % (m*n))/m];
92047c6ae99SBarry Smith         /* z_t = z; */
92147c6ae99SBarry Smith         s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
9228865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
92347c6ae99SBarry Smith       }
92447c6ae99SBarry Smith     }
92547c6ae99SBarry Smith   }
92647c6ae99SBarry Smith 
92747c6ae99SBarry Smith   /* Upper Level */
92847c6ae99SBarry Smith   for (k=0; k<s_z; k++) {
92947c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
93047c6ae99SBarry Smith       if (n18 >= 0) { /* left below */
931ce00eea3SSatish Balay         x_t = lx[n18 % m];
93247c6ae99SBarry Smith         y_t = ly[(n18 % (m*n))/m];
93347c6ae99SBarry Smith         /* z_t = lz[n18 / (m*n)]; */
93447c6ae99SBarry Smith         s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
9358865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n18] - (s_y-i)*x_t -s_x + x_t*y_t; /* 2d case */
9368865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
93747c6ae99SBarry Smith       }
93847c6ae99SBarry Smith       if (n19 >= 0) { /* directly below */
93947c6ae99SBarry Smith         x_t = x;
94047c6ae99SBarry Smith         y_t = ly[(n19 % (m*n))/m];
94147c6ae99SBarry Smith         /* z_t = lz[n19 / (m*n)]; */
94247c6ae99SBarry Smith         s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
9438865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n19] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */
9448865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
94547c6ae99SBarry Smith       }
94647c6ae99SBarry Smith       if (n20 >= 0) { /* right below */
947ce00eea3SSatish Balay         x_t = lx[n20 % m];
94847c6ae99SBarry Smith         y_t = ly[(n20 % (m*n))/m];
94947c6ae99SBarry Smith         /* z_t = lz[n20 / (m*n)]; */
95047c6ae99SBarry Smith         s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
9518865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n20] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */
9528865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
95347c6ae99SBarry Smith       }
95447c6ae99SBarry Smith     }
95547c6ae99SBarry Smith 
95647c6ae99SBarry Smith     for (i=0; i<y; i++) {
95747c6ae99SBarry Smith       if (n21 >= 0) { /* directly left */
958ce00eea3SSatish Balay         x_t = lx[n21 % m];
95947c6ae99SBarry Smith         y_t = y;
96047c6ae99SBarry Smith         /* z_t = lz[n21 / (m*n)]; */
96147c6ae99SBarry Smith         s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
9628865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n21] + (i+1)*x_t - s_x;  /* 2d case */
9638865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
96447c6ae99SBarry Smith       }
96547c6ae99SBarry Smith 
96647c6ae99SBarry Smith       if (n22 >= 0) { /* middle */
96747c6ae99SBarry Smith         x_t = x;
96847c6ae99SBarry Smith         y_t = y;
96947c6ae99SBarry Smith         /* z_t = lz[n22 / (m*n)]; */
97047c6ae99SBarry Smith         s_t = bases[n22] + i*x_t + k*x_t*y_t;
9718865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n22] + i*x_t; /* 2d case */
9728865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
973bff4a2f0SMatthew G. Knepley       } else if (bz == DM_BOUNDARY_MIRROR) {
9748865f1eaSKarl Rupp         for (j=0; j<x; j++) idx[nn++] = 0;
97547c6ae99SBarry Smith       }
97647c6ae99SBarry Smith 
97747c6ae99SBarry Smith       if (n23 >= 0) { /* directly right */
978ce00eea3SSatish Balay         x_t = lx[n23 % m];
97947c6ae99SBarry Smith         y_t = y;
98047c6ae99SBarry Smith         /* z_t = lz[n23 / (m*n)]; */
98147c6ae99SBarry Smith         s_t = bases[n23] + i*x_t + k*x_t*y_t;
9828865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n23] + i*x_t; /* 2d case */
9838865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
98447c6ae99SBarry Smith       }
98547c6ae99SBarry Smith     }
98647c6ae99SBarry Smith 
98747c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
98847c6ae99SBarry Smith       if (n24 >= 0) { /* left above */
989ce00eea3SSatish Balay         x_t = lx[n24 % m];
99047c6ae99SBarry Smith         y_t = ly[(n24 % (m*n))/m];
99147c6ae99SBarry Smith         /* z_t = lz[n24 / (m*n)]; */
99247c6ae99SBarry Smith         s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
9938865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n24] + i*x_t - s_x; /* 2d case */
9948865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
99547c6ae99SBarry Smith       }
99647c6ae99SBarry Smith       if (n25 >= 0) { /* directly above */
99747c6ae99SBarry Smith         x_t = x;
99847c6ae99SBarry Smith         y_t = ly[(n25 % (m*n))/m];
99947c6ae99SBarry Smith         /* z_t = lz[n25 / (m*n)]; */
100047c6ae99SBarry Smith         s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
10018865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n25] + (i-1)*x_t; /* 2d case */
10028865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
100347c6ae99SBarry Smith       }
100447c6ae99SBarry Smith       if (n26 >= 0) { /* right above */
1005ce00eea3SSatish Balay         x_t = lx[n26 % m];
100647c6ae99SBarry Smith         y_t = ly[(n26 % (m*n))/m];
100747c6ae99SBarry Smith         /* z_t = lz[n26 / (m*n)]; */
100847c6ae99SBarry Smith         s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
10098865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n26] + (i-1)*x_t; /* 2d case */
10108865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
101147c6ae99SBarry Smith       }
101247c6ae99SBarry Smith     }
101347c6ae99SBarry Smith   }
1014ce00eea3SSatish Balay 
1015b1fb7eb7SBarry Smith   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from);CHKERRQ(ierr);
101647c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
10173bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr);
1018fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
1019fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
102047c6ae99SBarry Smith 
1021d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
102247c6ae99SBarry Smith     n0  = sn0;  n1  = sn1;  n2  = sn2;  n3  = sn3;  n5  = sn5;  n6  = sn6; n7 = sn7;
102347c6ae99SBarry Smith     n8  = sn8;  n9  = sn9;  n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18;
102447c6ae99SBarry Smith     n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25;
102547c6ae99SBarry Smith     n26 = sn26;
1026ce00eea3SSatish Balay   }
102747c6ae99SBarry Smith 
102888661749SPeter Brune   if (((stencil_type == DMDA_STENCIL_STAR) ||
1029bff4a2f0SMatthew G. Knepley       (bx != DM_BOUNDARY_PERIODIC && bx) ||
1030bff4a2f0SMatthew G. Knepley       (by != DM_BOUNDARY_PERIODIC && by) ||
1031bff4a2f0SMatthew G. Knepley        (bz != DM_BOUNDARY_PERIODIC && bz))) {
1032ce00eea3SSatish Balay     /*
1033ce00eea3SSatish Balay         Recompute the local to global mappings, this time keeping the
1034ce00eea3SSatish Balay       information about the cross corner processor numbers.
1035ce00eea3SSatish Balay     */
103647c6ae99SBarry Smith     nn = 0;
103747c6ae99SBarry Smith     /* Bottom Level */
103847c6ae99SBarry Smith     for (k=0; k<s_z; k++) {
103947c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
104047c6ae99SBarry Smith         if (n0 >= 0) { /* left below */
1041ce00eea3SSatish Balay           x_t = lx[n0 % m];
104247c6ae99SBarry Smith           y_t = ly[(n0 % (m*n))/m];
104347c6ae99SBarry Smith           z_t = lz[n0 / (m*n)];
104447c6ae99SBarry 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;
10458865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1046ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Ys-ys < 0 && Zs-zs < 0) {
10478865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
104847c6ae99SBarry Smith         }
104947c6ae99SBarry Smith         if (n1 >= 0) { /* directly below */
105047c6ae99SBarry Smith           x_t = x;
105147c6ae99SBarry Smith           y_t = ly[(n1 % (m*n))/m];
105247c6ae99SBarry Smith           z_t = lz[n1 / (m*n)];
105347c6ae99SBarry 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;
10548865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1055ce00eea3SSatish Balay         } else if (Ys-ys < 0 && Zs-zs < 0) {
10568865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
105747c6ae99SBarry Smith         }
105847c6ae99SBarry Smith         if (n2 >= 0) { /* right below */
1059ce00eea3SSatish Balay           x_t = lx[n2 % m];
106047c6ae99SBarry Smith           y_t = ly[(n2 % (m*n))/m];
106147c6ae99SBarry Smith           z_t = lz[n2 / (m*n)];
106247c6ae99SBarry 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;
10638865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1064ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Ys-ys < 0 && Zs-zs < 0) {
10658865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
106647c6ae99SBarry Smith         }
106747c6ae99SBarry Smith       }
106847c6ae99SBarry Smith 
106947c6ae99SBarry Smith       for (i=0; i<y; i++) {
107047c6ae99SBarry Smith         if (n3 >= 0) { /* directly left */
1071ce00eea3SSatish Balay           x_t = lx[n3 % m];
107247c6ae99SBarry Smith           y_t = y;
107347c6ae99SBarry Smith           z_t = lz[n3 / (m*n)];
107447c6ae99SBarry Smith           s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
10758865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1076ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Zs-zs < 0) {
10778865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
107847c6ae99SBarry Smith         }
107947c6ae99SBarry Smith 
108047c6ae99SBarry Smith         if (n4 >= 0) { /* middle */
108147c6ae99SBarry Smith           x_t = x;
108247c6ae99SBarry Smith           y_t = y;
108347c6ae99SBarry Smith           z_t = lz[n4 / (m*n)];
108447c6ae99SBarry Smith           s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
10858865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1086ce00eea3SSatish Balay         } else if (Zs-zs < 0) {
1087bff4a2f0SMatthew G. Knepley           if (bz == DM_BOUNDARY_MIRROR) {
10888865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = 0;
1089c2859e5eSBarry Smith           } else {
10908865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = -1;
109147c6ae99SBarry Smith           }
1092c2859e5eSBarry Smith         }
109347c6ae99SBarry Smith 
109447c6ae99SBarry Smith         if (n5 >= 0) { /* directly right */
1095ce00eea3SSatish Balay           x_t = lx[n5 % m];
109647c6ae99SBarry Smith           y_t = y;
109747c6ae99SBarry Smith           z_t = lz[n5 / (m*n)];
109847c6ae99SBarry Smith           s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
10998865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1100ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Zs-zs < 0) {
11018865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
110247c6ae99SBarry Smith         }
110347c6ae99SBarry Smith       }
110447c6ae99SBarry Smith 
110547c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
110647c6ae99SBarry Smith         if (n6 >= 0) { /* left above */
1107ce00eea3SSatish Balay           x_t = lx[n6 % m];
110847c6ae99SBarry Smith           y_t = ly[(n6 % (m*n))/m];
110947c6ae99SBarry Smith           z_t = lz[n6 / (m*n)];
111047c6ae99SBarry Smith           s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
11118865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1112ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ye-Ye < 0 && Zs-zs < 0) {
11138865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
111447c6ae99SBarry Smith         }
111547c6ae99SBarry Smith         if (n7 >= 0) { /* directly above */
111647c6ae99SBarry Smith           x_t = x;
111747c6ae99SBarry Smith           y_t = ly[(n7 % (m*n))/m];
111847c6ae99SBarry Smith           z_t = lz[n7 / (m*n)];
111947c6ae99SBarry Smith           s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
11208865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1121ce00eea3SSatish Balay         } else if (ye-Ye < 0 && Zs-zs < 0) {
11228865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
112347c6ae99SBarry Smith         }
112447c6ae99SBarry Smith         if (n8 >= 0) { /* right above */
1125ce00eea3SSatish Balay           x_t = lx[n8 % m];
112647c6ae99SBarry Smith           y_t = ly[(n8 % (m*n))/m];
112747c6ae99SBarry Smith           z_t = lz[n8 / (m*n)];
112847c6ae99SBarry Smith           s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
11298865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1130ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ye-Ye < 0 && Zs-zs < 0) {
11318865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
113247c6ae99SBarry Smith         }
113347c6ae99SBarry Smith       }
113447c6ae99SBarry Smith     }
113547c6ae99SBarry Smith 
113647c6ae99SBarry Smith     /* Middle Level */
113747c6ae99SBarry Smith     for (k=0; k<z; k++) {
113847c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
113947c6ae99SBarry Smith         if (n9 >= 0) { /* left below */
1140ce00eea3SSatish Balay           x_t = lx[n9 % m];
114147c6ae99SBarry Smith           y_t = ly[(n9 % (m*n))/m];
114247c6ae99SBarry Smith           /* z_t = z; */
114347c6ae99SBarry Smith           s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
11448865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1145ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Ys-ys < 0) {
11468865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
114747c6ae99SBarry Smith         }
114847c6ae99SBarry Smith         if (n10 >= 0) { /* directly below */
114947c6ae99SBarry Smith           x_t = x;
115047c6ae99SBarry Smith           y_t = ly[(n10 % (m*n))/m];
115147c6ae99SBarry Smith           /* z_t = z; */
115247c6ae99SBarry Smith           s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
11538865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1154ce00eea3SSatish Balay         } else if (Ys-ys < 0) {
1155bff4a2f0SMatthew G. Knepley           if (by == DM_BOUNDARY_MIRROR) {
11568865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = -1;
1157c2859e5eSBarry Smith           } else {
11588865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = -1;
1159c2859e5eSBarry Smith           }
116047c6ae99SBarry Smith         }
116147c6ae99SBarry Smith         if (n11 >= 0) { /* right below */
1162ce00eea3SSatish Balay           x_t = lx[n11 % m];
116347c6ae99SBarry Smith           y_t = ly[(n11 % (m*n))/m];
116447c6ae99SBarry Smith           /* z_t = z; */
116547c6ae99SBarry Smith           s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
11668865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1167ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Ys-ys < 0) {
11688865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
116947c6ae99SBarry Smith         }
117047c6ae99SBarry Smith       }
117147c6ae99SBarry Smith 
117247c6ae99SBarry Smith       for (i=0; i<y; i++) {
117347c6ae99SBarry Smith         if (n12 >= 0) { /* directly left */
1174ce00eea3SSatish Balay           x_t = lx[n12 % m];
117547c6ae99SBarry Smith           y_t = y;
117647c6ae99SBarry Smith           /* z_t = z; */
117747c6ae99SBarry Smith           s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
11788865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1179ce00eea3SSatish Balay         } else if (Xs-xs < 0) {
1180bff4a2f0SMatthew G. Knepley           if (bx == DM_BOUNDARY_MIRROR) {
11818865f1eaSKarl Rupp             for (j=0; j<s_x; j++) idx[nn++] = 0;
1182c2859e5eSBarry Smith           } else {
11838865f1eaSKarl Rupp             for (j=0; j<s_x; j++) idx[nn++] = -1;
118447c6ae99SBarry Smith           }
1185c2859e5eSBarry Smith         }
118647c6ae99SBarry Smith 
118747c6ae99SBarry Smith         /* Interior */
118847c6ae99SBarry Smith         s_t = bases[rank] + i*x + k*x*y;
11898865f1eaSKarl Rupp         for (j=0; j<x; j++) idx[nn++] = s_t++;
119047c6ae99SBarry Smith 
119147c6ae99SBarry Smith         if (n14 >= 0) { /* directly right */
1192ce00eea3SSatish Balay           x_t = lx[n14 % m];
119347c6ae99SBarry Smith           y_t = y;
119447c6ae99SBarry Smith           /* z_t = z; */
119547c6ae99SBarry Smith           s_t = bases[n14] + i*x_t + k*x_t*y_t;
11968865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1197ce00eea3SSatish Balay         } else if (xe-Xe < 0) {
1198bff4a2f0SMatthew G. Knepley           if (bx == DM_BOUNDARY_MIRROR) {
11998865f1eaSKarl Rupp             for (j=0; j<s_x; j++) idx[nn++] = 0;
1200c2859e5eSBarry Smith           } else {
12018865f1eaSKarl Rupp             for (j=0; j<s_x; j++) idx[nn++] = -1;
120247c6ae99SBarry Smith           }
120347c6ae99SBarry Smith         }
1204c2859e5eSBarry Smith       }
120547c6ae99SBarry Smith 
120647c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
120747c6ae99SBarry Smith         if (n15 >= 0) { /* left above */
1208ce00eea3SSatish Balay           x_t = lx[n15 % m];
120947c6ae99SBarry Smith           y_t = ly[(n15 % (m*n))/m];
121047c6ae99SBarry Smith           /* z_t = z; */
121147c6ae99SBarry Smith           s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
12128865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1213ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ye-Ye < 0) {
12148865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
121547c6ae99SBarry Smith         }
121647c6ae99SBarry Smith         if (n16 >= 0) { /* directly above */
121747c6ae99SBarry Smith           x_t = x;
121847c6ae99SBarry Smith           y_t = ly[(n16 % (m*n))/m];
121947c6ae99SBarry Smith           /* z_t = z; */
122047c6ae99SBarry Smith           s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
12218865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1222ce00eea3SSatish Balay         } else if (ye-Ye < 0) {
1223bff4a2f0SMatthew G. Knepley           if (by == DM_BOUNDARY_MIRROR) {
12248865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = 0;
1225c2859e5eSBarry Smith           } else {
12268865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = -1;
122747c6ae99SBarry Smith           }
1228c2859e5eSBarry Smith         }
122947c6ae99SBarry Smith         if (n17 >= 0) { /* right above */
1230ce00eea3SSatish Balay           x_t = lx[n17 % m];
123147c6ae99SBarry Smith           y_t = ly[(n17 % (m*n))/m];
123247c6ae99SBarry Smith           /* z_t = z; */
123347c6ae99SBarry Smith           s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
12348865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1235ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ye-Ye < 0) {
12368865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
123747c6ae99SBarry Smith         }
123847c6ae99SBarry Smith       }
123947c6ae99SBarry Smith     }
124047c6ae99SBarry Smith 
124147c6ae99SBarry Smith     /* Upper Level */
124247c6ae99SBarry Smith     for (k=0; k<s_z; k++) {
124347c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
124447c6ae99SBarry Smith         if (n18 >= 0) { /* left below */
1245ce00eea3SSatish Balay           x_t = lx[n18 % m];
124647c6ae99SBarry Smith           y_t = ly[(n18 % (m*n))/m];
124747c6ae99SBarry Smith           /* z_t = lz[n18 / (m*n)]; */
124847c6ae99SBarry Smith           s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
12498865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1250ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Ys-ys < 0 && ze-Ze < 0) {
12518865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
125247c6ae99SBarry Smith         }
125347c6ae99SBarry Smith         if (n19 >= 0) { /* directly below */
125447c6ae99SBarry Smith           x_t = x;
125547c6ae99SBarry Smith           y_t = ly[(n19 % (m*n))/m];
125647c6ae99SBarry Smith           /* z_t = lz[n19 / (m*n)]; */
125747c6ae99SBarry Smith           s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
12588865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1259ce00eea3SSatish Balay         } else if (Ys-ys < 0 && ze-Ze < 0) {
12608865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
126147c6ae99SBarry Smith         }
126247c6ae99SBarry Smith         if (n20 >= 0) { /* right below */
1263ce00eea3SSatish Balay           x_t = lx[n20 % m];
126447c6ae99SBarry Smith           y_t = ly[(n20 % (m*n))/m];
126547c6ae99SBarry Smith           /* z_t = lz[n20 / (m*n)]; */
126647c6ae99SBarry Smith           s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
12678865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1268ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Ys-ys < 0 && ze-Ze < 0) {
12698865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
127047c6ae99SBarry Smith         }
127147c6ae99SBarry Smith       }
127247c6ae99SBarry Smith 
127347c6ae99SBarry Smith       for (i=0; i<y; i++) {
127447c6ae99SBarry Smith         if (n21 >= 0) { /* directly left */
1275ce00eea3SSatish Balay           x_t = lx[n21 % m];
127647c6ae99SBarry Smith           y_t = y;
127747c6ae99SBarry Smith           /* z_t = lz[n21 / (m*n)]; */
127847c6ae99SBarry Smith           s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
12798865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1280ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ze-Ze < 0) {
12818865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
128247c6ae99SBarry Smith         }
128347c6ae99SBarry Smith 
128447c6ae99SBarry Smith         if (n22 >= 0) { /* middle */
128547c6ae99SBarry Smith           x_t = x;
128647c6ae99SBarry Smith           y_t = y;
128747c6ae99SBarry Smith           /* z_t = lz[n22 / (m*n)]; */
128847c6ae99SBarry Smith           s_t = bases[n22] + i*x_t + k*x_t*y_t;
12898865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1290ce00eea3SSatish Balay         } else if (ze-Ze < 0) {
1291bff4a2f0SMatthew G. Knepley           if (bz == DM_BOUNDARY_MIRROR) {
12928865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = 0;
1293c2859e5eSBarry Smith           } else {
12948865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = -1;
129547c6ae99SBarry Smith           }
1296c2859e5eSBarry Smith         }
129747c6ae99SBarry Smith 
129847c6ae99SBarry Smith         if (n23 >= 0) { /* directly right */
1299ce00eea3SSatish Balay           x_t = lx[n23 % m];
130047c6ae99SBarry Smith           y_t = y;
130147c6ae99SBarry Smith           /* z_t = lz[n23 / (m*n)]; */
130247c6ae99SBarry Smith           s_t = bases[n23] + i*x_t + k*x_t*y_t;
13038865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1304ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ze-Ze < 0) {
13058865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
130647c6ae99SBarry Smith         }
130747c6ae99SBarry Smith       }
130847c6ae99SBarry Smith 
130947c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
131047c6ae99SBarry Smith         if (n24 >= 0) { /* left above */
1311ce00eea3SSatish Balay           x_t = lx[n24 % m];
131247c6ae99SBarry Smith           y_t = ly[(n24 % (m*n))/m];
131347c6ae99SBarry Smith           /* z_t = lz[n24 / (m*n)]; */
131447c6ae99SBarry Smith           s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
13158865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1316ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ye-Ye < 0 && ze-Ze < 0) {
13178865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
131847c6ae99SBarry Smith         }
131947c6ae99SBarry Smith         if (n25 >= 0) { /* directly above */
132047c6ae99SBarry Smith           x_t = x;
132147c6ae99SBarry Smith           y_t = ly[(n25 % (m*n))/m];
132247c6ae99SBarry Smith           /* z_t = lz[n25 / (m*n)]; */
132347c6ae99SBarry Smith           s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
13248865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1325ce00eea3SSatish Balay         } else if (ye-Ye < 0 && ze-Ze < 0) {
13268865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
132747c6ae99SBarry Smith         }
132847c6ae99SBarry Smith         if (n26 >= 0) { /* right above */
1329ce00eea3SSatish Balay           x_t = lx[n26 % m];
133047c6ae99SBarry Smith           y_t = ly[(n26 % (m*n))/m];
133147c6ae99SBarry Smith           /* z_t = lz[n26 / (m*n)]; */
133247c6ae99SBarry Smith           s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
13338865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1334ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ye-Ye < 0 && ze-Ze < 0) {
13358865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
133647c6ae99SBarry Smith         }
133747c6ae99SBarry Smith       }
133847c6ae99SBarry Smith     }
133947c6ae99SBarry Smith   }
134047c6ae99SBarry Smith   /*
134147c6ae99SBarry Smith      Set the local to global ordering in the global vector, this allows use
134247c6ae99SBarry Smith      of VecSetValuesLocal().
134347c6ae99SBarry Smith   */
134445b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr);
13453bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr);
134647c6ae99SBarry Smith 
1347ce00eea3SSatish Balay   ierr  = PetscFree2(bases,ldims);CHKERRQ(ierr);
1348ce00eea3SSatish Balay   dd->m = m;  dd->n  = n;  dd->p  = p;
1349ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
1350ce00eea3SSatish Balay   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = zs; dd->ze = ze;
1351ce00eea3SSatish Balay   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = Zs; dd->Ze = Ze;
1352ce00eea3SSatish Balay 
1353fcfd50ebSBarry Smith   ierr = VecDestroy(&local);CHKERRQ(ierr);
1354fcfd50ebSBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
1355ce00eea3SSatish Balay 
1356ce00eea3SSatish Balay   dd->gtol      = gtol;
1357ce00eea3SSatish Balay   dd->base      = base;
1358ce00eea3SSatish Balay   da->ops->view = DMView_DA_3d;
13590298fd71SBarry Smith   dd->ltol      = NULL;
13600298fd71SBarry Smith   dd->ao        = NULL;
136147c6ae99SBarry Smith   PetscFunctionReturn(0);
136247c6ae99SBarry Smith }
1363564755cdSBarry Smith 
136447c6ae99SBarry Smith 
136547c6ae99SBarry Smith #undef __FUNCT__
1366aa219208SBarry Smith #define __FUNCT__ "DMDACreate3d"
136747c6ae99SBarry Smith /*@C
1368aa219208SBarry Smith    DMDACreate3d - Creates an object that will manage the communication of three-dimensional
136947c6ae99SBarry Smith    regular array data that is distributed across some processors.
137047c6ae99SBarry Smith 
137147c6ae99SBarry Smith    Collective on MPI_Comm
137247c6ae99SBarry Smith 
137347c6ae99SBarry Smith    Input Parameters:
137447c6ae99SBarry Smith +  comm - MPI communicator
13751321219cSEthan Coon .  bx,by,bz - type of ghost nodes the array have.
1376bff4a2f0SMatthew G. Knepley          Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC.
1377aa219208SBarry Smith .  stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX)
137847c6ae99SBarry 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
137947c6ae99SBarry Smith             from the command line with -da_grid_x <M> -da_grid_y <N> -da_grid_z <P>)
138047c6ae99SBarry Smith .  m,n,p - corresponding number of processors in each dimension
138147c6ae99SBarry Smith            (or PETSC_DECIDE to have calculated)
138247c6ae99SBarry Smith .  dof - number of degrees of freedom per node
138310d7c030SMatthew G Knepley .  s - stencil width
138410d7c030SMatthew G Knepley -  lx, ly, lz - arrays containing the number of nodes in each cell along
13850298fd71SBarry Smith           the x, y, and z coordinates, or NULL. If non-null, these
138647c6ae99SBarry Smith           must be of length as m,n,p and the corresponding
138747c6ae99SBarry Smith           m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of
138847c6ae99SBarry Smith           the ly[] must N, sum of the lz[] must be P
138947c6ae99SBarry Smith 
139047c6ae99SBarry Smith    Output Parameter:
139147c6ae99SBarry Smith .  da - the resulting distributed array object
139247c6ae99SBarry Smith 
139347c6ae99SBarry Smith    Options Database Key:
1394706a11cbSBarry Smith +  -dm_view - Calls DMView() at the conclusion of DMDACreate3d()
139547c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
139647c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
139747c6ae99SBarry Smith .  -da_grid_z <nz> - number of grid points in z direction, if P < 0
139847c6ae99SBarry Smith .  -da_processors_x <MX> - number of processors in x direction
139947c6ae99SBarry Smith .  -da_processors_y <MY> - number of processors in y direction
140047c6ae99SBarry Smith .  -da_processors_z <MZ> - number of processors in z direction
1401e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement ratio in x direction
1402e0f5d30fSBarry Smith .  -da_refine_y <ry> - refinement ratio in y direction
1403e0f5d30fSBarry Smith .  -da_refine_z <rz>- refinement ratio in z directio
1404e0f5d30fSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating it, , if M, N, or P < 0
140547c6ae99SBarry Smith 
140647c6ae99SBarry Smith    Level: beginner
140747c6ae99SBarry Smith 
140847c6ae99SBarry Smith    Notes:
1409aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
1410aa219208SBarry Smith    standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
141147c6ae99SBarry Smith    the standard 27-pt stencil.
141247c6ae99SBarry Smith 
1413aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
1414564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
1415564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
141647c6ae99SBarry Smith 
141747c6ae99SBarry Smith .keywords: distributed array, create, three-dimensional
141847c6ae99SBarry Smith 
1419aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate2d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
142099f0b487SRichard Tran Mills           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(),
1421d461ba97SBarry Smith           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
142247c6ae99SBarry Smith 
142347c6ae99SBarry Smith @*/
1424bff4a2f0SMatthew G. Knepley PetscErrorCode  DMDACreate3d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMBoundaryType bz,DMDAStencilType stencil_type,PetscInt M,
14259a42bb27SBarry 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)
142647c6ae99SBarry Smith {
142747c6ae99SBarry Smith   PetscErrorCode ierr;
142847c6ae99SBarry Smith 
142947c6ae99SBarry Smith   PetscFunctionBegin;
1430aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
1431c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*da, 3);CHKERRQ(ierr);
1432aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, P);CHKERRQ(ierr);
1433aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, p);CHKERRQ(ierr);
14341321219cSEthan Coon   ierr = DMDASetBoundaryType(*da, bx, by, bz);CHKERRQ(ierr);
1435aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
1436aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
1437aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
1438aa219208SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, lz);CHKERRQ(ierr);
143947c6ae99SBarry Smith   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
14409a42bb27SBarry Smith   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
14419a42bb27SBarry Smith   ierr = DMSetUp(*da);CHKERRQ(ierr);
144247c6ae99SBarry Smith   PetscFunctionReturn(0);
144347c6ae99SBarry Smith }
1444