xref: /petsc/src/dm/impls/da/da3.c (revision 45b6f7e9ecd564e8bb535880b270d4a8084ff211)
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 
74035e84dSBarry 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"
129a42bb27SBarry Smith 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 
347b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);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);
547b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);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 */
1298ea3bf28SBarry Smith         base = (dd->Xe-dd->Xs)*(dd->Ye-dd->Ys)*(k-dd->Zs);
130*45b6f7e9SBarry Smith         ierr = ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr);
13147c6ae99SBarry Smith         plane=k;
13247c6ae99SBarry Smith         /* Keep z wrap around points on the dradrawg */
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);
15147c6ae99SBarry Smith             base+=dd->w;
15247c6ae99SBarry Smith           }
15347c6ae99SBarry Smith         }
154*45b6f7e9SBarry 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;
192ce00eea3SSatish Balay   PetscInt         Xs,Xe,Ys,Ye,Zs,Ze,IXs,IXe,IYs,IYe,IZs,IZe,start,end,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;
20147c6ae99SBarry Smith   VecScatter       ltog,gtol;
202*45b6f7e9SBarry 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;
425778a2246SBarry Smith   ierr       = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr);
426ce00eea3SSatish Balay   dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs)*dof;
427778a2246SBarry Smith   ierr       = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,0,&local);CHKERRQ(ierr);
42847c6ae99SBarry Smith 
42947c6ae99SBarry Smith   /* generate appropriate vector scatters */
43047c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
43147c6ae99SBarry Smith   ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr);
432ce00eea3SSatish Balay   ierr = ISCreateStride(comm,x*y*z*dof,start,1,&to);CHKERRQ(ierr);
43347c6ae99SBarry Smith 
434785e854fSJed Brown   ierr   = PetscMalloc1(x*y*z,&idx);CHKERRQ(ierr);
435ce00eea3SSatish Balay   left   = xs - Xs; right = left + x;
43647c6ae99SBarry Smith   bottom = ys - Ys; top = bottom + y;
43747c6ae99SBarry Smith   down   = zs - Zs; up  = down + z;
43847c6ae99SBarry Smith   count  = 0;
43947c6ae99SBarry Smith   for (i=down; i<up; i++) {
44047c6ae99SBarry Smith     for (j=bottom; j<top; j++) {
441ce00eea3SSatish Balay       for (k=left; k<right; k++) {
442ce00eea3SSatish Balay         idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
44347c6ae99SBarry Smith       }
44447c6ae99SBarry Smith     }
44547c6ae99SBarry Smith   }
44647c6ae99SBarry Smith 
44747c6ae99SBarry Smith   ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
44847c6ae99SBarry Smith   ierr = VecScatterCreate(local,from,global,to,&ltog);CHKERRQ(ierr);
4493bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)ltog);CHKERRQ(ierr);
450fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
451fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
45247c6ae99SBarry Smith 
453ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
454ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
455d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_BOX) {
456db87c5ecSEthan Coon     count = (IXe-IXs)*(IYe-IYs)*(IZe-IZs);
457785e854fSJed Brown     ierr  = PetscMalloc1(count,&idx);CHKERRQ(ierr);
458ce00eea3SSatish Balay 
459ce00eea3SSatish Balay     left   = IXs - Xs; right = left + (IXe-IXs);
460ce00eea3SSatish Balay     bottom = IYs - Ys; top = bottom + (IYe-IYs);
461ce00eea3SSatish Balay     down   = IZs - Zs; up  = down + (IZe-IZs);
462ce00eea3SSatish Balay     count  = 0;
463ce00eea3SSatish Balay     for (i=down; i<up; i++) {
464ce00eea3SSatish Balay       for (j=bottom; j<top; j++) {
465ce00eea3SSatish Balay         for (k=left; k<right; k++) {
466ce00eea3SSatish Balay           idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
467ce00eea3SSatish Balay         }
468ce00eea3SSatish Balay       }
469ce00eea3SSatish Balay     }
470ce00eea3SSatish Balay     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
471ce00eea3SSatish Balay 
47247c6ae99SBarry Smith   } else {
47347c6ae99SBarry Smith     /* This is way ugly! We need to list the funny cross type region */
474db87c5ecSEthan Coon     count = ((ys-IYs) + (IYe-ye))*x*z + ((xs-IXs) + (IXe-xe))*y*z + ((zs-IZs) + (IZe-ze))*x*y + x*y*z;
475785e854fSJed Brown     ierr  = PetscMalloc1(count,&idx);CHKERRQ(ierr);
476ce00eea3SSatish Balay 
477ce00eea3SSatish Balay     left   = xs - Xs; right = left + x;
47847c6ae99SBarry Smith     bottom = ys - Ys; top = bottom + y;
47947c6ae99SBarry Smith     down   = zs - Zs;   up  = down + z;
48047c6ae99SBarry Smith     count  = 0;
481ce00eea3SSatish Balay     /* the bottom chunck */
482ce00eea3SSatish Balay     for (i=(IZs-Zs); i<down; i++) {
48347c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
484ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
48547c6ae99SBarry Smith       }
48647c6ae99SBarry Smith     }
48747c6ae99SBarry Smith     /* the middle piece */
48847c6ae99SBarry Smith     for (i=down; i<up; i++) {
48947c6ae99SBarry Smith       /* front */
490ce00eea3SSatish Balay       for (j=(IYs-Ys); j<bottom; j++) {
491ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
49247c6ae99SBarry Smith       }
49347c6ae99SBarry Smith       /* middle */
49447c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
495ce00eea3SSatish Balay         for (k=IXs-Xs; k<IXe-Xs; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
49647c6ae99SBarry Smith       }
49747c6ae99SBarry Smith       /* back */
498ce00eea3SSatish Balay       for (j=top; j<top+IYe-ye; j++) {
499ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
50047c6ae99SBarry Smith       }
50147c6ae99SBarry Smith     }
50247c6ae99SBarry Smith     /* the top piece */
503ce00eea3SSatish Balay     for (i=up; i<up+IZe-ze; i++) {
50447c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
505ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
50647c6ae99SBarry Smith       }
50747c6ae99SBarry Smith     }
50847c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
50947c6ae99SBarry Smith   }
51047c6ae99SBarry Smith 
51147c6ae99SBarry Smith   /* determine who lies on each side of use stored in    n24 n25 n26
51247c6ae99SBarry Smith                                                          n21 n22 n23
51347c6ae99SBarry Smith                                                          n18 n19 n20
51447c6ae99SBarry Smith 
51547c6ae99SBarry Smith                                                          n15 n16 n17
51647c6ae99SBarry Smith                                                          n12     n14
51747c6ae99SBarry Smith                                                          n9  n10 n11
51847c6ae99SBarry Smith 
51947c6ae99SBarry Smith                                                          n6  n7  n8
52047c6ae99SBarry Smith                                                          n3  n4  n5
52147c6ae99SBarry Smith                                                          n0  n1  n2
52247c6ae99SBarry Smith   */
52347c6ae99SBarry Smith 
52447c6ae99SBarry Smith   /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
52547c6ae99SBarry Smith   /* Assume Nodes are Internal to the Cube */
52647c6ae99SBarry Smith   n0 = rank - m*n - m - 1;
52747c6ae99SBarry Smith   n1 = rank - m*n - m;
52847c6ae99SBarry Smith   n2 = rank - m*n - m + 1;
52947c6ae99SBarry Smith   n3 = rank - m*n -1;
53047c6ae99SBarry Smith   n4 = rank - m*n;
53147c6ae99SBarry Smith   n5 = rank - m*n + 1;
53247c6ae99SBarry Smith   n6 = rank - m*n + m - 1;
53347c6ae99SBarry Smith   n7 = rank - m*n + m;
53447c6ae99SBarry Smith   n8 = rank - m*n + m + 1;
53547c6ae99SBarry Smith 
53647c6ae99SBarry Smith   n9  = rank - m - 1;
53747c6ae99SBarry Smith   n10 = rank - m;
53847c6ae99SBarry Smith   n11 = rank - m + 1;
53947c6ae99SBarry Smith   n12 = rank - 1;
54047c6ae99SBarry Smith   n14 = rank + 1;
54147c6ae99SBarry Smith   n15 = rank + m - 1;
54247c6ae99SBarry Smith   n16 = rank + m;
54347c6ae99SBarry Smith   n17 = rank + m + 1;
54447c6ae99SBarry Smith 
54547c6ae99SBarry Smith   n18 = rank + m*n - m - 1;
54647c6ae99SBarry Smith   n19 = rank + m*n - m;
54747c6ae99SBarry Smith   n20 = rank + m*n - m + 1;
54847c6ae99SBarry Smith   n21 = rank + m*n - 1;
54947c6ae99SBarry Smith   n22 = rank + m*n;
55047c6ae99SBarry Smith   n23 = rank + m*n + 1;
55147c6ae99SBarry Smith   n24 = rank + m*n + m - 1;
55247c6ae99SBarry Smith   n25 = rank + m*n + m;
55347c6ae99SBarry Smith   n26 = rank + m*n + m + 1;
55447c6ae99SBarry Smith 
55547c6ae99SBarry Smith   /* Assume Pieces are on Faces of Cube */
55647c6ae99SBarry Smith 
55747c6ae99SBarry Smith   if (xs == 0) { /* First assume not corner or edge */
55847c6ae99SBarry Smith     n0  = rank       -1 - (m*n);
55947c6ae99SBarry Smith     n3  = rank + m   -1 - (m*n);
56047c6ae99SBarry Smith     n6  = rank + 2*m -1 - (m*n);
56147c6ae99SBarry Smith     n9  = rank       -1;
56247c6ae99SBarry Smith     n12 = rank + m   -1;
56347c6ae99SBarry Smith     n15 = rank + 2*m -1;
56447c6ae99SBarry Smith     n18 = rank       -1 + (m*n);
56547c6ae99SBarry Smith     n21 = rank + m   -1 + (m*n);
56647c6ae99SBarry Smith     n24 = rank + 2*m -1 + (m*n);
56747c6ae99SBarry Smith   }
56847c6ae99SBarry Smith 
569ce00eea3SSatish Balay   if (xe == M) { /* First assume not corner or edge */
57047c6ae99SBarry Smith     n2  = rank -2*m +1 - (m*n);
57147c6ae99SBarry Smith     n5  = rank - m  +1 - (m*n);
57247c6ae99SBarry Smith     n8  = rank      +1 - (m*n);
57347c6ae99SBarry Smith     n11 = rank -2*m +1;
57447c6ae99SBarry Smith     n14 = rank - m  +1;
57547c6ae99SBarry Smith     n17 = rank      +1;
57647c6ae99SBarry Smith     n20 = rank -2*m +1 + (m*n);
57747c6ae99SBarry Smith     n23 = rank - m  +1 + (m*n);
57847c6ae99SBarry Smith     n26 = rank      +1 + (m*n);
57947c6ae99SBarry Smith   }
58047c6ae99SBarry Smith 
58147c6ae99SBarry Smith   if (ys==0) { /* First assume not corner or edge */
58247c6ae99SBarry Smith     n0  = rank + m * (n-1) -1 - (m*n);
58347c6ae99SBarry Smith     n1  = rank + m * (n-1)    - (m*n);
58447c6ae99SBarry Smith     n2  = rank + m * (n-1) +1 - (m*n);
58547c6ae99SBarry Smith     n9  = rank + m * (n-1) -1;
58647c6ae99SBarry Smith     n10 = rank + m * (n-1);
58747c6ae99SBarry Smith     n11 = rank + m * (n-1) +1;
58847c6ae99SBarry Smith     n18 = rank + m * (n-1) -1 + (m*n);
58947c6ae99SBarry Smith     n19 = rank + m * (n-1)    + (m*n);
59047c6ae99SBarry Smith     n20 = rank + m * (n-1) +1 + (m*n);
59147c6ae99SBarry Smith   }
59247c6ae99SBarry Smith 
59347c6ae99SBarry Smith   if (ye == N) { /* First assume not corner or edge */
59447c6ae99SBarry Smith     n6  = rank - m * (n-1) -1 - (m*n);
59547c6ae99SBarry Smith     n7  = rank - m * (n-1)    - (m*n);
59647c6ae99SBarry Smith     n8  = rank - m * (n-1) +1 - (m*n);
59747c6ae99SBarry Smith     n15 = rank - m * (n-1) -1;
59847c6ae99SBarry Smith     n16 = rank - m * (n-1);
59947c6ae99SBarry Smith     n17 = rank - m * (n-1) +1;
60047c6ae99SBarry Smith     n24 = rank - m * (n-1) -1 + (m*n);
60147c6ae99SBarry Smith     n25 = rank - m * (n-1)    + (m*n);
60247c6ae99SBarry Smith     n26 = rank - m * (n-1) +1 + (m*n);
60347c6ae99SBarry Smith   }
60447c6ae99SBarry Smith 
60547c6ae99SBarry Smith   if (zs == 0) { /* First assume not corner or edge */
60647c6ae99SBarry Smith     n0 = size - (m*n) + rank - m - 1;
60747c6ae99SBarry Smith     n1 = size - (m*n) + rank - m;
60847c6ae99SBarry Smith     n2 = size - (m*n) + rank - m + 1;
60947c6ae99SBarry Smith     n3 = size - (m*n) + rank - 1;
61047c6ae99SBarry Smith     n4 = size - (m*n) + rank;
61147c6ae99SBarry Smith     n5 = size - (m*n) + rank + 1;
61247c6ae99SBarry Smith     n6 = size - (m*n) + rank + m - 1;
61347c6ae99SBarry Smith     n7 = size - (m*n) + rank + m;
61447c6ae99SBarry Smith     n8 = size - (m*n) + rank + m + 1;
61547c6ae99SBarry Smith   }
61647c6ae99SBarry Smith 
61747c6ae99SBarry Smith   if (ze == P) { /* First assume not corner or edge */
61847c6ae99SBarry Smith     n18 = (m*n) - (size-rank) - m - 1;
61947c6ae99SBarry Smith     n19 = (m*n) - (size-rank) - m;
62047c6ae99SBarry Smith     n20 = (m*n) - (size-rank) - m + 1;
62147c6ae99SBarry Smith     n21 = (m*n) - (size-rank) - 1;
62247c6ae99SBarry Smith     n22 = (m*n) - (size-rank);
62347c6ae99SBarry Smith     n23 = (m*n) - (size-rank) + 1;
62447c6ae99SBarry Smith     n24 = (m*n) - (size-rank) + m - 1;
62547c6ae99SBarry Smith     n25 = (m*n) - (size-rank) + m;
62647c6ae99SBarry Smith     n26 = (m*n) - (size-rank) + m + 1;
62747c6ae99SBarry Smith   }
62847c6ae99SBarry Smith 
62947c6ae99SBarry Smith   if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */
63047c6ae99SBarry Smith     n0 = size - m*n + rank + m-1 - m;
63147c6ae99SBarry Smith     n3 = size - m*n + rank + m-1;
63247c6ae99SBarry Smith     n6 = size - m*n + rank + m-1 + m;
63347c6ae99SBarry Smith   }
63447c6ae99SBarry Smith 
63547c6ae99SBarry Smith   if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */
63647c6ae99SBarry Smith     n18 = m*n - (size - rank) + m-1 - m;
63747c6ae99SBarry Smith     n21 = m*n - (size - rank) + m-1;
63847c6ae99SBarry Smith     n24 = m*n - (size - rank) + m-1 + m;
63947c6ae99SBarry Smith   }
64047c6ae99SBarry Smith 
64147c6ae99SBarry Smith   if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */
64247c6ae99SBarry Smith     n0  = rank + m*n -1 - m*n;
64347c6ae99SBarry Smith     n9  = rank + m*n -1;
64447c6ae99SBarry Smith     n18 = rank + m*n -1 + m*n;
64547c6ae99SBarry Smith   }
64647c6ae99SBarry Smith 
64747c6ae99SBarry Smith   if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */
64847c6ae99SBarry Smith     n6  = rank - m*(n-1) + m-1 - m*n;
64947c6ae99SBarry Smith     n15 = rank - m*(n-1) + m-1;
65047c6ae99SBarry Smith     n24 = rank - m*(n-1) + m-1 + m*n;
65147c6ae99SBarry Smith   }
65247c6ae99SBarry Smith 
653ce00eea3SSatish Balay   if ((xe==M) && (zs==0)) { /* Assume an edge, not corner */
65447c6ae99SBarry Smith     n2 = size - (m*n-rank) - (m-1) - m;
65547c6ae99SBarry Smith     n5 = size - (m*n-rank) - (m-1);
65647c6ae99SBarry Smith     n8 = size - (m*n-rank) - (m-1) + m;
65747c6ae99SBarry Smith   }
65847c6ae99SBarry Smith 
659ce00eea3SSatish Balay   if ((xe==M) && (ze==P)) { /* Assume an edge, not corner */
66047c6ae99SBarry Smith     n20 = m*n - (size - rank) - (m-1) - m;
66147c6ae99SBarry Smith     n23 = m*n - (size - rank) - (m-1);
66247c6ae99SBarry Smith     n26 = m*n - (size - rank) - (m-1) + m;
66347c6ae99SBarry Smith   }
66447c6ae99SBarry Smith 
665ce00eea3SSatish Balay   if ((xe==M) && (ys==0)) { /* Assume an edge, not corner */
66647c6ae99SBarry Smith     n2  = rank + m*(n-1) - (m-1) - m*n;
66747c6ae99SBarry Smith     n11 = rank + m*(n-1) - (m-1);
66847c6ae99SBarry Smith     n20 = rank + m*(n-1) - (m-1) + m*n;
66947c6ae99SBarry Smith   }
67047c6ae99SBarry Smith 
671ce00eea3SSatish Balay   if ((xe==M) && (ye==N)) { /* Assume an edge, not corner */
67247c6ae99SBarry Smith     n8  = rank - m*n +1 - m*n;
67347c6ae99SBarry Smith     n17 = rank - m*n +1;
67447c6ae99SBarry Smith     n26 = rank - m*n +1 + m*n;
67547c6ae99SBarry Smith   }
67647c6ae99SBarry Smith 
67747c6ae99SBarry Smith   if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */
67847c6ae99SBarry Smith     n0 = size - m + rank -1;
67947c6ae99SBarry Smith     n1 = size - m + rank;
68047c6ae99SBarry Smith     n2 = size - m + rank +1;
68147c6ae99SBarry Smith   }
68247c6ae99SBarry Smith 
68347c6ae99SBarry Smith   if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */
68447c6ae99SBarry Smith     n18 = m*n - (size - rank) + m*(n-1) -1;
68547c6ae99SBarry Smith     n19 = m*n - (size - rank) + m*(n-1);
68647c6ae99SBarry Smith     n20 = m*n - (size - rank) + m*(n-1) +1;
68747c6ae99SBarry Smith   }
68847c6ae99SBarry Smith 
68947c6ae99SBarry Smith   if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */
69047c6ae99SBarry Smith     n6 = size - (m*n-rank) - m * (n-1) -1;
69147c6ae99SBarry Smith     n7 = size - (m*n-rank) - m * (n-1);
69247c6ae99SBarry Smith     n8 = size - (m*n-rank) - m * (n-1) +1;
69347c6ae99SBarry Smith   }
69447c6ae99SBarry Smith 
69547c6ae99SBarry Smith   if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */
69647c6ae99SBarry Smith     n24 = rank - (size-m) -1;
69747c6ae99SBarry Smith     n25 = rank - (size-m);
69847c6ae99SBarry Smith     n26 = rank - (size-m) +1;
69947c6ae99SBarry Smith   }
70047c6ae99SBarry Smith 
70147c6ae99SBarry Smith   /* Check for Corners */
7028865f1eaSKarl Rupp   if ((xs==0) && (ys==0) && (zs==0)) n0  = size -1;
7038865f1eaSKarl Rupp   if ((xs==0) && (ys==0) && (ze==P)) n18 = m*n-1;
7048865f1eaSKarl Rupp   if ((xs==0) && (ye==N) && (zs==0)) n6  = (size-1)-m*(n-1);
7058865f1eaSKarl Rupp   if ((xs==0) && (ye==N) && (ze==P)) n24 = m-1;
7068865f1eaSKarl Rupp   if ((xe==M) && (ys==0) && (zs==0)) n2  = size-m;
7078865f1eaSKarl Rupp   if ((xe==M) && (ys==0) && (ze==P)) n20 = m*n-m;
7088865f1eaSKarl Rupp   if ((xe==M) && (ye==N) && (zs==0)) n8  = size-m*n;
7098865f1eaSKarl Rupp   if ((xe==M) && (ye==N) && (ze==P)) n26 = 0;
71047c6ae99SBarry Smith 
71147c6ae99SBarry Smith   /* Check for when not X,Y, and Z Periodic */
71247c6ae99SBarry Smith 
71347c6ae99SBarry Smith   /* If not X periodic */
714bff4a2f0SMatthew G. Knepley   if (bx != DM_BOUNDARY_PERIODIC) {
7158865f1eaSKarl Rupp     if (xs==0) n0 = n3 = n6 = n9  = n12 = n15 = n18 = n21 = n24 = -2;
7168865f1eaSKarl Rupp     if (xe==M) n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2;
71747c6ae99SBarry Smith   }
71847c6ae99SBarry Smith 
71947c6ae99SBarry Smith   /* If not Y periodic */
720bff4a2f0SMatthew G. Knepley   if (by != DM_BOUNDARY_PERIODIC) {
7218865f1eaSKarl Rupp     if (ys==0) n0 = n1 = n2 = n9  = n10 = n11 = n18 = n19 = n20 = -2;
7228865f1eaSKarl Rupp     if (ye==N) n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2;
72347c6ae99SBarry Smith   }
72447c6ae99SBarry Smith 
72547c6ae99SBarry Smith   /* If not Z periodic */
726bff4a2f0SMatthew G. Knepley   if (bz != DM_BOUNDARY_PERIODIC) {
7278865f1eaSKarl Rupp     if (zs==0) n0  = n1  = n2  = n3  = n4  = n5  = n6  = n7  = n8  = -2;
7288865f1eaSKarl Rupp     if (ze==P) n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;
72947c6ae99SBarry Smith   }
73047c6ae99SBarry Smith 
731785e854fSJed Brown   ierr = PetscMalloc1(27,&dd->neighbors);CHKERRQ(ierr);
7328865f1eaSKarl Rupp 
73347c6ae99SBarry Smith   dd->neighbors[0]  = n0;
73447c6ae99SBarry Smith   dd->neighbors[1]  = n1;
73547c6ae99SBarry Smith   dd->neighbors[2]  = n2;
73647c6ae99SBarry Smith   dd->neighbors[3]  = n3;
73747c6ae99SBarry Smith   dd->neighbors[4]  = n4;
73847c6ae99SBarry Smith   dd->neighbors[5]  = n5;
73947c6ae99SBarry Smith   dd->neighbors[6]  = n6;
74047c6ae99SBarry Smith   dd->neighbors[7]  = n7;
74147c6ae99SBarry Smith   dd->neighbors[8]  = n8;
74247c6ae99SBarry Smith   dd->neighbors[9]  = n9;
74347c6ae99SBarry Smith   dd->neighbors[10] = n10;
74447c6ae99SBarry Smith   dd->neighbors[11] = n11;
74547c6ae99SBarry Smith   dd->neighbors[12] = n12;
74647c6ae99SBarry Smith   dd->neighbors[13] = rank;
74747c6ae99SBarry Smith   dd->neighbors[14] = n14;
74847c6ae99SBarry Smith   dd->neighbors[15] = n15;
74947c6ae99SBarry Smith   dd->neighbors[16] = n16;
75047c6ae99SBarry Smith   dd->neighbors[17] = n17;
75147c6ae99SBarry Smith   dd->neighbors[18] = n18;
75247c6ae99SBarry Smith   dd->neighbors[19] = n19;
75347c6ae99SBarry Smith   dd->neighbors[20] = n20;
75447c6ae99SBarry Smith   dd->neighbors[21] = n21;
75547c6ae99SBarry Smith   dd->neighbors[22] = n22;
75647c6ae99SBarry Smith   dd->neighbors[23] = n23;
75747c6ae99SBarry Smith   dd->neighbors[24] = n24;
75847c6ae99SBarry Smith   dd->neighbors[25] = n25;
75947c6ae99SBarry Smith   dd->neighbors[26] = n26;
76047c6ae99SBarry Smith 
76147c6ae99SBarry Smith   /* If star stencil then delete the corner neighbors */
762d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
76347c6ae99SBarry Smith     /* save information about corner neighbors */
76447c6ae99SBarry Smith     sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7;
76547c6ae99SBarry Smith     sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18;
76647c6ae99SBarry Smith     sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25;
76747c6ae99SBarry Smith     sn26 = n26;
7688865f1eaSKarl Rupp     n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 = n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1;
76947c6ae99SBarry Smith   }
77047c6ae99SBarry Smith 
771785e854fSJed Brown   ierr = PetscMalloc1((Xe-Xs)*(Ye-Ys)*(Ze-Zs),&idx);CHKERRQ(ierr);
7723bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)da,(Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt));CHKERRQ(ierr);
77347c6ae99SBarry Smith 
77447c6ae99SBarry Smith   nn = 0;
77547c6ae99SBarry Smith   /* Bottom Level */
77647c6ae99SBarry Smith   for (k=0; k<s_z; k++) {
77747c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
77847c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
779ce00eea3SSatish Balay         x_t = lx[n0 % m];
78047c6ae99SBarry Smith         y_t = ly[(n0 % (m*n))/m];
78147c6ae99SBarry Smith         z_t = lz[n0 / (m*n)];
78247c6ae99SBarry 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;
7838865f1eaSKarl 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 */
7848865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
78547c6ae99SBarry Smith       }
78647c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
78747c6ae99SBarry Smith         x_t = x;
78847c6ae99SBarry Smith         y_t = ly[(n1 % (m*n))/m];
78947c6ae99SBarry Smith         z_t = lz[n1 / (m*n)];
79047c6ae99SBarry 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;
7918865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */
7928865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
79347c6ae99SBarry Smith       }
79447c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
795ce00eea3SSatish Balay         x_t = lx[n2 % m];
79647c6ae99SBarry Smith         y_t = ly[(n2 % (m*n))/m];
79747c6ae99SBarry Smith         z_t = lz[n2 / (m*n)];
79847c6ae99SBarry 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;
7998865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t; /* 2D case */
8008865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
80147c6ae99SBarry Smith       }
80247c6ae99SBarry Smith     }
80347c6ae99SBarry Smith 
80447c6ae99SBarry Smith     for (i=0; i<y; i++) {
80547c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
806ce00eea3SSatish Balay         x_t = lx[n3 % m];
80747c6ae99SBarry Smith         y_t = y;
80847c6ae99SBarry Smith         z_t = lz[n3 / (m*n)];
80947c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8108865f1eaSKarl 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 */
8118865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
81247c6ae99SBarry Smith       }
81347c6ae99SBarry Smith 
81447c6ae99SBarry Smith       if (n4 >= 0) { /* middle */
81547c6ae99SBarry Smith         x_t = x;
81647c6ae99SBarry Smith         y_t = y;
81747c6ae99SBarry Smith         z_t = lz[n4 / (m*n)];
81847c6ae99SBarry Smith         s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8198865f1eaSKarl 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 */
8208865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
821bff4a2f0SMatthew G. Knepley       } else if (bz == DM_BOUNDARY_MIRROR) {
8228865f1eaSKarl Rupp         for (j=0; j<x; j++) idx[nn++] = 0;
82347c6ae99SBarry Smith       }
82447c6ae99SBarry Smith 
82547c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
826ce00eea3SSatish Balay         x_t = lx[n5 % m];
82747c6ae99SBarry Smith         y_t = y;
82847c6ae99SBarry Smith         z_t = lz[n5 / (m*n)];
82947c6ae99SBarry Smith         s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8308865f1eaSKarl 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 */
8318865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
83247c6ae99SBarry Smith       }
83347c6ae99SBarry Smith     }
83447c6ae99SBarry Smith 
83547c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
83647c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
837ce00eea3SSatish Balay         x_t = lx[n6 % m];
83847c6ae99SBarry Smith         y_t = ly[(n6 % (m*n))/m];
83947c6ae99SBarry Smith         z_t = lz[n6 / (m*n)];
84047c6ae99SBarry Smith         s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8418865f1eaSKarl 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 */
8428865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
84347c6ae99SBarry Smith       }
84447c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
84547c6ae99SBarry Smith         x_t = x;
84647c6ae99SBarry Smith         y_t = ly[(n7 % (m*n))/m];
84747c6ae99SBarry Smith         z_t = lz[n7 / (m*n)];
84847c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8498865f1eaSKarl 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 */
8508865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
85147c6ae99SBarry Smith       }
85247c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
853ce00eea3SSatish Balay         x_t = lx[n8 % m];
85447c6ae99SBarry Smith         y_t = ly[(n8 % (m*n))/m];
85547c6ae99SBarry Smith         z_t = lz[n8 / (m*n)];
85647c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
8578865f1eaSKarl 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 */
8588865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
85947c6ae99SBarry Smith       }
86047c6ae99SBarry Smith     }
86147c6ae99SBarry Smith   }
86247c6ae99SBarry Smith 
86347c6ae99SBarry Smith   /* Middle Level */
86447c6ae99SBarry Smith   for (k=0; k<z; k++) {
86547c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
86647c6ae99SBarry Smith       if (n9 >= 0) { /* left below */
867ce00eea3SSatish Balay         x_t = lx[n9 % m];
86847c6ae99SBarry Smith         y_t = ly[(n9 % (m*n))/m];
86947c6ae99SBarry Smith         /* z_t = z; */
87047c6ae99SBarry Smith         s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
8718865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
87247c6ae99SBarry Smith       }
87347c6ae99SBarry Smith       if (n10 >= 0) { /* directly below */
87447c6ae99SBarry Smith         x_t = x;
87547c6ae99SBarry Smith         y_t = ly[(n10 % (m*n))/m];
87647c6ae99SBarry Smith         /* z_t = z; */
87747c6ae99SBarry Smith         s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
8788865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
879bff4a2f0SMatthew G. Knepley       }  else if (by == DM_BOUNDARY_MIRROR) {
8808865f1eaSKarl Rupp         for (j=0; j<x; j++) idx[nn++] = 0;
88147c6ae99SBarry Smith       }
88247c6ae99SBarry Smith       if (n11 >= 0) { /* right below */
883ce00eea3SSatish Balay         x_t = lx[n11 % m];
88447c6ae99SBarry Smith         y_t = ly[(n11 % (m*n))/m];
88547c6ae99SBarry Smith         /* z_t = z; */
88647c6ae99SBarry Smith         s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
8878865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
88847c6ae99SBarry Smith       }
88947c6ae99SBarry Smith     }
89047c6ae99SBarry Smith 
89147c6ae99SBarry Smith     for (i=0; i<y; i++) {
89247c6ae99SBarry Smith       if (n12 >= 0) { /* directly left */
893ce00eea3SSatish Balay         x_t = lx[n12 % m];
89447c6ae99SBarry Smith         y_t = y;
89547c6ae99SBarry Smith         /* z_t = z; */
89647c6ae99SBarry Smith         s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
8978865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
898bff4a2f0SMatthew G. Knepley       }  else if (bx == DM_BOUNDARY_MIRROR) {
8998865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = 0;
90047c6ae99SBarry Smith       }
90147c6ae99SBarry Smith 
90247c6ae99SBarry Smith       /* Interior */
90347c6ae99SBarry Smith       s_t = bases[rank] + i*x + k*x*y;
9048865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = s_t++;
90547c6ae99SBarry Smith 
90647c6ae99SBarry Smith       if (n14 >= 0) { /* directly right */
907ce00eea3SSatish Balay         x_t = lx[n14 % m];
90847c6ae99SBarry Smith         y_t = y;
90947c6ae99SBarry Smith         /* z_t = z; */
91047c6ae99SBarry Smith         s_t = bases[n14] + i*x_t + k*x_t*y_t;
9118865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
912bff4a2f0SMatthew G. Knepley       } else if (bx == DM_BOUNDARY_MIRROR) {
9138865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = 0;
91447c6ae99SBarry Smith       }
91547c6ae99SBarry Smith     }
91647c6ae99SBarry Smith 
91747c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
91847c6ae99SBarry Smith       if (n15 >= 0) { /* left above */
919ce00eea3SSatish Balay         x_t = lx[n15 % m];
92047c6ae99SBarry Smith         y_t = ly[(n15 % (m*n))/m];
92147c6ae99SBarry Smith         /* z_t = z; */
92247c6ae99SBarry Smith         s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
9238865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
92447c6ae99SBarry Smith       }
92547c6ae99SBarry Smith       if (n16 >= 0) { /* directly above */
92647c6ae99SBarry Smith         x_t = x;
92747c6ae99SBarry Smith         y_t = ly[(n16 % (m*n))/m];
92847c6ae99SBarry Smith         /* z_t = z; */
92947c6ae99SBarry Smith         s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
9308865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
931bff4a2f0SMatthew G. Knepley       } else if (by == DM_BOUNDARY_MIRROR) {
9328865f1eaSKarl Rupp         for (j=0; j<x; j++) idx[nn++] = 0;
93347c6ae99SBarry Smith       }
93447c6ae99SBarry Smith       if (n17 >= 0) { /* right above */
935ce00eea3SSatish Balay         x_t = lx[n17 % m];
93647c6ae99SBarry Smith         y_t = ly[(n17 % (m*n))/m];
93747c6ae99SBarry Smith         /* z_t = z; */
93847c6ae99SBarry Smith         s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
9398865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
94047c6ae99SBarry Smith       }
94147c6ae99SBarry Smith     }
94247c6ae99SBarry Smith   }
94347c6ae99SBarry Smith 
94447c6ae99SBarry Smith   /* Upper Level */
94547c6ae99SBarry Smith   for (k=0; k<s_z; k++) {
94647c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
94747c6ae99SBarry Smith       if (n18 >= 0) { /* left below */
948ce00eea3SSatish Balay         x_t = lx[n18 % m];
94947c6ae99SBarry Smith         y_t = ly[(n18 % (m*n))/m];
95047c6ae99SBarry Smith         /* z_t = lz[n18 / (m*n)]; */
95147c6ae99SBarry Smith         s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
9528865f1eaSKarl 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 */
9538865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
95447c6ae99SBarry Smith       }
95547c6ae99SBarry Smith       if (n19 >= 0) { /* directly below */
95647c6ae99SBarry Smith         x_t = x;
95747c6ae99SBarry Smith         y_t = ly[(n19 % (m*n))/m];
95847c6ae99SBarry Smith         /* z_t = lz[n19 / (m*n)]; */
95947c6ae99SBarry Smith         s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
9608865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n19] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */
9618865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
96247c6ae99SBarry Smith       }
96347c6ae99SBarry Smith       if (n20 >= 0) { /* right below */
964ce00eea3SSatish Balay         x_t = lx[n20 % m];
96547c6ae99SBarry Smith         y_t = ly[(n20 % (m*n))/m];
96647c6ae99SBarry Smith         /* z_t = lz[n20 / (m*n)]; */
96747c6ae99SBarry Smith         s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
9688865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n20] - (s_y+1-i)*x_t + x_t*y_t; /* 2d case */
9698865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
97047c6ae99SBarry Smith       }
97147c6ae99SBarry Smith     }
97247c6ae99SBarry Smith 
97347c6ae99SBarry Smith     for (i=0; i<y; i++) {
97447c6ae99SBarry Smith       if (n21 >= 0) { /* directly left */
975ce00eea3SSatish Balay         x_t = lx[n21 % m];
97647c6ae99SBarry Smith         y_t = y;
97747c6ae99SBarry Smith         /* z_t = lz[n21 / (m*n)]; */
97847c6ae99SBarry Smith         s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
9798865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n21] + (i+1)*x_t - s_x;  /* 2d case */
9808865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
98147c6ae99SBarry Smith       }
98247c6ae99SBarry Smith 
98347c6ae99SBarry Smith       if (n22 >= 0) { /* middle */
98447c6ae99SBarry Smith         x_t = x;
98547c6ae99SBarry Smith         y_t = y;
98647c6ae99SBarry Smith         /* z_t = lz[n22 / (m*n)]; */
98747c6ae99SBarry Smith         s_t = bases[n22] + i*x_t + k*x_t*y_t;
9888865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n22] + i*x_t; /* 2d case */
9898865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
990bff4a2f0SMatthew G. Knepley       } else if (bz == DM_BOUNDARY_MIRROR) {
9918865f1eaSKarl Rupp         for (j=0; j<x; j++) idx[nn++] = 0;
99247c6ae99SBarry Smith       }
99347c6ae99SBarry Smith 
99447c6ae99SBarry Smith       if (n23 >= 0) { /* directly right */
995ce00eea3SSatish Balay         x_t = lx[n23 % m];
99647c6ae99SBarry Smith         y_t = y;
99747c6ae99SBarry Smith         /* z_t = lz[n23 / (m*n)]; */
99847c6ae99SBarry Smith         s_t = bases[n23] + i*x_t + k*x_t*y_t;
9998865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n23] + i*x_t; /* 2d case */
10008865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
100147c6ae99SBarry Smith       }
100247c6ae99SBarry Smith     }
100347c6ae99SBarry Smith 
100447c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
100547c6ae99SBarry Smith       if (n24 >= 0) { /* left above */
1006ce00eea3SSatish Balay         x_t = lx[n24 % m];
100747c6ae99SBarry Smith         y_t = ly[(n24 % (m*n))/m];
100847c6ae99SBarry Smith         /* z_t = lz[n24 / (m*n)]; */
100947c6ae99SBarry Smith         s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
10108865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n24] + i*x_t - s_x; /* 2d case */
10118865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
101247c6ae99SBarry Smith       }
101347c6ae99SBarry Smith       if (n25 >= 0) { /* directly above */
101447c6ae99SBarry Smith         x_t = x;
101547c6ae99SBarry Smith         y_t = ly[(n25 % (m*n))/m];
101647c6ae99SBarry Smith         /* z_t = lz[n25 / (m*n)]; */
101747c6ae99SBarry Smith         s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
10188865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n25] + (i-1)*x_t; /* 2d case */
10198865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
102047c6ae99SBarry Smith       }
102147c6ae99SBarry Smith       if (n26 >= 0) { /* right above */
1022ce00eea3SSatish Balay         x_t = lx[n26 % m];
102347c6ae99SBarry Smith         y_t = ly[(n26 % (m*n))/m];
102447c6ae99SBarry Smith         /* z_t = lz[n26 / (m*n)]; */
102547c6ae99SBarry Smith         s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
10268865f1eaSKarl Rupp         if (twod && (s_t >= M*N*P)) s_t = bases[n26] + (i-1)*x_t; /* 2d case */
10278865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
102847c6ae99SBarry Smith       }
102947c6ae99SBarry Smith     }
103047c6ae99SBarry Smith   }
1031ce00eea3SSatish Balay 
1032ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
103347c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
10343bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr);
1035fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
1036fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
103747c6ae99SBarry Smith 
1038d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
103947c6ae99SBarry Smith     n0  = sn0;  n1  = sn1;  n2  = sn2;  n3  = sn3;  n5  = sn5;  n6  = sn6; n7 = sn7;
104047c6ae99SBarry Smith     n8  = sn8;  n9  = sn9;  n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18;
104147c6ae99SBarry Smith     n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25;
104247c6ae99SBarry Smith     n26 = sn26;
1043ce00eea3SSatish Balay   }
104447c6ae99SBarry Smith 
104588661749SPeter Brune   if (((stencil_type == DMDA_STENCIL_STAR) ||
1046bff4a2f0SMatthew G. Knepley       (bx != DM_BOUNDARY_PERIODIC && bx) ||
1047bff4a2f0SMatthew G. Knepley       (by != DM_BOUNDARY_PERIODIC && by) ||
1048bff4a2f0SMatthew G. Knepley        (bz != DM_BOUNDARY_PERIODIC && bz))) {
1049ce00eea3SSatish Balay     /*
1050ce00eea3SSatish Balay         Recompute the local to global mappings, this time keeping the
1051ce00eea3SSatish Balay       information about the cross corner processor numbers.
1052ce00eea3SSatish Balay     */
105347c6ae99SBarry Smith     nn = 0;
105447c6ae99SBarry Smith     /* Bottom Level */
105547c6ae99SBarry Smith     for (k=0; k<s_z; k++) {
105647c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
105747c6ae99SBarry Smith         if (n0 >= 0) { /* left below */
1058ce00eea3SSatish Balay           x_t = lx[n0 % m];
105947c6ae99SBarry Smith           y_t = ly[(n0 % (m*n))/m];
106047c6ae99SBarry Smith           z_t = lz[n0 / (m*n)];
106147c6ae99SBarry 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;
10628865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1063ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Ys-ys < 0 && Zs-zs < 0) {
10648865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
106547c6ae99SBarry Smith         }
106647c6ae99SBarry Smith         if (n1 >= 0) { /* directly below */
106747c6ae99SBarry Smith           x_t = x;
106847c6ae99SBarry Smith           y_t = ly[(n1 % (m*n))/m];
106947c6ae99SBarry Smith           z_t = lz[n1 / (m*n)];
107047c6ae99SBarry 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;
10718865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1072ce00eea3SSatish Balay         } else if (Ys-ys < 0 && Zs-zs < 0) {
10738865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
107447c6ae99SBarry Smith         }
107547c6ae99SBarry Smith         if (n2 >= 0) { /* right below */
1076ce00eea3SSatish Balay           x_t = lx[n2 % m];
107747c6ae99SBarry Smith           y_t = ly[(n2 % (m*n))/m];
107847c6ae99SBarry Smith           z_t = lz[n2 / (m*n)];
107947c6ae99SBarry 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;
10808865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1081ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Ys-ys < 0 && Zs-zs < 0) {
10828865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
108347c6ae99SBarry Smith         }
108447c6ae99SBarry Smith       }
108547c6ae99SBarry Smith 
108647c6ae99SBarry Smith       for (i=0; i<y; i++) {
108747c6ae99SBarry Smith         if (n3 >= 0) { /* directly left */
1088ce00eea3SSatish Balay           x_t = lx[n3 % m];
108947c6ae99SBarry Smith           y_t = y;
109047c6ae99SBarry Smith           z_t = lz[n3 / (m*n)];
109147c6ae99SBarry Smith           s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
10928865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1093ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Zs-zs < 0) {
10948865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
109547c6ae99SBarry Smith         }
109647c6ae99SBarry Smith 
109747c6ae99SBarry Smith         if (n4 >= 0) { /* middle */
109847c6ae99SBarry Smith           x_t = x;
109947c6ae99SBarry Smith           y_t = y;
110047c6ae99SBarry Smith           z_t = lz[n4 / (m*n)];
110147c6ae99SBarry Smith           s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
11028865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1103ce00eea3SSatish Balay         } else if (Zs-zs < 0) {
1104bff4a2f0SMatthew G. Knepley           if (bz == DM_BOUNDARY_MIRROR) {
11058865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = 0;
1106c2859e5eSBarry Smith           } else {
11078865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = -1;
110847c6ae99SBarry Smith           }
1109c2859e5eSBarry Smith         }
111047c6ae99SBarry Smith 
111147c6ae99SBarry Smith         if (n5 >= 0) { /* directly right */
1112ce00eea3SSatish Balay           x_t = lx[n5 % m];
111347c6ae99SBarry Smith           y_t = y;
111447c6ae99SBarry Smith           z_t = lz[n5 / (m*n)];
111547c6ae99SBarry Smith           s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
11168865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1117ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Zs-zs < 0) {
11188865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
111947c6ae99SBarry Smith         }
112047c6ae99SBarry Smith       }
112147c6ae99SBarry Smith 
112247c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
112347c6ae99SBarry Smith         if (n6 >= 0) { /* left above */
1124ce00eea3SSatish Balay           x_t = lx[n6 % m];
112547c6ae99SBarry Smith           y_t = ly[(n6 % (m*n))/m];
112647c6ae99SBarry Smith           z_t = lz[n6 / (m*n)];
112747c6ae99SBarry Smith           s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
11288865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1129ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ye-Ye < 0 && Zs-zs < 0) {
11308865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
113147c6ae99SBarry Smith         }
113247c6ae99SBarry Smith         if (n7 >= 0) { /* directly above */
113347c6ae99SBarry Smith           x_t = x;
113447c6ae99SBarry Smith           y_t = ly[(n7 % (m*n))/m];
113547c6ae99SBarry Smith           z_t = lz[n7 / (m*n)];
113647c6ae99SBarry Smith           s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
11378865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1138ce00eea3SSatish Balay         } else if (ye-Ye < 0 && Zs-zs < 0) {
11398865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
114047c6ae99SBarry Smith         }
114147c6ae99SBarry Smith         if (n8 >= 0) { /* right above */
1142ce00eea3SSatish Balay           x_t = lx[n8 % m];
114347c6ae99SBarry Smith           y_t = ly[(n8 % (m*n))/m];
114447c6ae99SBarry Smith           z_t = lz[n8 / (m*n)];
114547c6ae99SBarry Smith           s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
11468865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1147ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ye-Ye < 0 && Zs-zs < 0) {
11488865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
114947c6ae99SBarry Smith         }
115047c6ae99SBarry Smith       }
115147c6ae99SBarry Smith     }
115247c6ae99SBarry Smith 
115347c6ae99SBarry Smith     /* Middle Level */
115447c6ae99SBarry Smith     for (k=0; k<z; k++) {
115547c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
115647c6ae99SBarry Smith         if (n9 >= 0) { /* left below */
1157ce00eea3SSatish Balay           x_t = lx[n9 % m];
115847c6ae99SBarry Smith           y_t = ly[(n9 % (m*n))/m];
115947c6ae99SBarry Smith           /* z_t = z; */
116047c6ae99SBarry Smith           s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
11618865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1162ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Ys-ys < 0) {
11638865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
116447c6ae99SBarry Smith         }
116547c6ae99SBarry Smith         if (n10 >= 0) { /* directly below */
116647c6ae99SBarry Smith           x_t = x;
116747c6ae99SBarry Smith           y_t = ly[(n10 % (m*n))/m];
116847c6ae99SBarry Smith           /* z_t = z; */
116947c6ae99SBarry Smith           s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
11708865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1171ce00eea3SSatish Balay         } else if (Ys-ys < 0) {
1172bff4a2f0SMatthew G. Knepley           if (by == DM_BOUNDARY_MIRROR) {
11738865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = -1;
1174c2859e5eSBarry Smith           } else {
11758865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = -1;
1176c2859e5eSBarry Smith           }
117747c6ae99SBarry Smith         }
117847c6ae99SBarry Smith         if (n11 >= 0) { /* right below */
1179ce00eea3SSatish Balay           x_t = lx[n11 % m];
118047c6ae99SBarry Smith           y_t = ly[(n11 % (m*n))/m];
118147c6ae99SBarry Smith           /* z_t = z; */
118247c6ae99SBarry Smith           s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
11838865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1184ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Ys-ys < 0) {
11858865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
118647c6ae99SBarry Smith         }
118747c6ae99SBarry Smith       }
118847c6ae99SBarry Smith 
118947c6ae99SBarry Smith       for (i=0; i<y; i++) {
119047c6ae99SBarry Smith         if (n12 >= 0) { /* directly left */
1191ce00eea3SSatish Balay           x_t = lx[n12 % m];
119247c6ae99SBarry Smith           y_t = y;
119347c6ae99SBarry Smith           /* z_t = z; */
119447c6ae99SBarry Smith           s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
11958865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1196ce00eea3SSatish Balay         } else if (Xs-xs < 0) {
1197bff4a2f0SMatthew G. Knepley           if (bx == DM_BOUNDARY_MIRROR) {
11988865f1eaSKarl Rupp             for (j=0; j<s_x; j++) idx[nn++] = 0;
1199c2859e5eSBarry Smith           } else {
12008865f1eaSKarl Rupp             for (j=0; j<s_x; j++) idx[nn++] = -1;
120147c6ae99SBarry Smith           }
1202c2859e5eSBarry Smith         }
120347c6ae99SBarry Smith 
120447c6ae99SBarry Smith         /* Interior */
120547c6ae99SBarry Smith         s_t = bases[rank] + i*x + k*x*y;
12068865f1eaSKarl Rupp         for (j=0; j<x; j++) idx[nn++] = s_t++;
120747c6ae99SBarry Smith 
120847c6ae99SBarry Smith         if (n14 >= 0) { /* directly right */
1209ce00eea3SSatish Balay           x_t = lx[n14 % m];
121047c6ae99SBarry Smith           y_t = y;
121147c6ae99SBarry Smith           /* z_t = z; */
121247c6ae99SBarry Smith           s_t = bases[n14] + i*x_t + k*x_t*y_t;
12138865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1214ce00eea3SSatish Balay         } else if (xe-Xe < 0) {
1215bff4a2f0SMatthew G. Knepley           if (bx == DM_BOUNDARY_MIRROR) {
12168865f1eaSKarl Rupp             for (j=0; j<s_x; j++) idx[nn++] = 0;
1217c2859e5eSBarry Smith           } else {
12188865f1eaSKarl Rupp             for (j=0; j<s_x; j++) idx[nn++] = -1;
121947c6ae99SBarry Smith           }
122047c6ae99SBarry Smith         }
1221c2859e5eSBarry Smith       }
122247c6ae99SBarry Smith 
122347c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
122447c6ae99SBarry Smith         if (n15 >= 0) { /* left above */
1225ce00eea3SSatish Balay           x_t = lx[n15 % m];
122647c6ae99SBarry Smith           y_t = ly[(n15 % (m*n))/m];
122747c6ae99SBarry Smith           /* z_t = z; */
122847c6ae99SBarry Smith           s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
12298865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1230ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ye-Ye < 0) {
12318865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
123247c6ae99SBarry Smith         }
123347c6ae99SBarry Smith         if (n16 >= 0) { /* directly above */
123447c6ae99SBarry Smith           x_t = x;
123547c6ae99SBarry Smith           y_t = ly[(n16 % (m*n))/m];
123647c6ae99SBarry Smith           /* z_t = z; */
123747c6ae99SBarry Smith           s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
12388865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1239ce00eea3SSatish Balay         } else if (ye-Ye < 0) {
1240bff4a2f0SMatthew G. Knepley           if (by == DM_BOUNDARY_MIRROR) {
12418865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = 0;
1242c2859e5eSBarry Smith           } else {
12438865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = -1;
124447c6ae99SBarry Smith           }
1245c2859e5eSBarry Smith         }
124647c6ae99SBarry Smith         if (n17 >= 0) { /* right above */
1247ce00eea3SSatish Balay           x_t = lx[n17 % m];
124847c6ae99SBarry Smith           y_t = ly[(n17 % (m*n))/m];
124947c6ae99SBarry Smith           /* z_t = z; */
125047c6ae99SBarry Smith           s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
12518865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1252ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ye-Ye < 0) {
12538865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
125447c6ae99SBarry Smith         }
125547c6ae99SBarry Smith       }
125647c6ae99SBarry Smith     }
125747c6ae99SBarry Smith 
125847c6ae99SBarry Smith     /* Upper Level */
125947c6ae99SBarry Smith     for (k=0; k<s_z; k++) {
126047c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
126147c6ae99SBarry Smith         if (n18 >= 0) { /* left below */
1262ce00eea3SSatish Balay           x_t = lx[n18 % m];
126347c6ae99SBarry Smith           y_t = ly[(n18 % (m*n))/m];
126447c6ae99SBarry Smith           /* z_t = lz[n18 / (m*n)]; */
126547c6ae99SBarry Smith           s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
12668865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1267ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Ys-ys < 0 && ze-Ze < 0) {
12688865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
126947c6ae99SBarry Smith         }
127047c6ae99SBarry Smith         if (n19 >= 0) { /* directly below */
127147c6ae99SBarry Smith           x_t = x;
127247c6ae99SBarry Smith           y_t = ly[(n19 % (m*n))/m];
127347c6ae99SBarry Smith           /* z_t = lz[n19 / (m*n)]; */
127447c6ae99SBarry Smith           s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
12758865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1276ce00eea3SSatish Balay         } else if (Ys-ys < 0 && ze-Ze < 0) {
12778865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
127847c6ae99SBarry Smith         }
127947c6ae99SBarry Smith         if (n20 >= 0) { /* right below */
1280ce00eea3SSatish Balay           x_t = lx[n20 % m];
128147c6ae99SBarry Smith           y_t = ly[(n20 % (m*n))/m];
128247c6ae99SBarry Smith           /* z_t = lz[n20 / (m*n)]; */
128347c6ae99SBarry Smith           s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
12848865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1285ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Ys-ys < 0 && ze-Ze < 0) {
12868865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
128747c6ae99SBarry Smith         }
128847c6ae99SBarry Smith       }
128947c6ae99SBarry Smith 
129047c6ae99SBarry Smith       for (i=0; i<y; i++) {
129147c6ae99SBarry Smith         if (n21 >= 0) { /* directly left */
1292ce00eea3SSatish Balay           x_t = lx[n21 % m];
129347c6ae99SBarry Smith           y_t = y;
129447c6ae99SBarry Smith           /* z_t = lz[n21 / (m*n)]; */
129547c6ae99SBarry Smith           s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
12968865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1297ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ze-Ze < 0) {
12988865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
129947c6ae99SBarry Smith         }
130047c6ae99SBarry Smith 
130147c6ae99SBarry Smith         if (n22 >= 0) { /* middle */
130247c6ae99SBarry Smith           x_t = x;
130347c6ae99SBarry Smith           y_t = y;
130447c6ae99SBarry Smith           /* z_t = lz[n22 / (m*n)]; */
130547c6ae99SBarry Smith           s_t = bases[n22] + i*x_t + k*x_t*y_t;
13068865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1307ce00eea3SSatish Balay         } else if (ze-Ze < 0) {
1308bff4a2f0SMatthew G. Knepley           if (bz == DM_BOUNDARY_MIRROR) {
13098865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = 0;
1310c2859e5eSBarry Smith           } else {
13118865f1eaSKarl Rupp             for (j=0; j<x; j++) idx[nn++] = -1;
131247c6ae99SBarry Smith           }
1313c2859e5eSBarry Smith         }
131447c6ae99SBarry Smith 
131547c6ae99SBarry Smith         if (n23 >= 0) { /* directly right */
1316ce00eea3SSatish Balay           x_t = lx[n23 % m];
131747c6ae99SBarry Smith           y_t = y;
131847c6ae99SBarry Smith           /* z_t = lz[n23 / (m*n)]; */
131947c6ae99SBarry Smith           s_t = bases[n23] + i*x_t + k*x_t*y_t;
13208865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1321ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ze-Ze < 0) {
13228865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
132347c6ae99SBarry Smith         }
132447c6ae99SBarry Smith       }
132547c6ae99SBarry Smith 
132647c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
132747c6ae99SBarry Smith         if (n24 >= 0) { /* left above */
1328ce00eea3SSatish Balay           x_t = lx[n24 % m];
132947c6ae99SBarry Smith           y_t = ly[(n24 % (m*n))/m];
133047c6ae99SBarry Smith           /* z_t = lz[n24 / (m*n)]; */
133147c6ae99SBarry Smith           s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
13328865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1333ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ye-Ye < 0 && ze-Ze < 0) {
13348865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
133547c6ae99SBarry Smith         }
133647c6ae99SBarry Smith         if (n25 >= 0) { /* directly above */
133747c6ae99SBarry Smith           x_t = x;
133847c6ae99SBarry Smith           y_t = ly[(n25 % (m*n))/m];
133947c6ae99SBarry Smith           /* z_t = lz[n25 / (m*n)]; */
134047c6ae99SBarry Smith           s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
13418865f1eaSKarl Rupp           for (j=0; j<x_t; j++) idx[nn++] = s_t++;
1342ce00eea3SSatish Balay         } else if (ye-Ye < 0 && ze-Ze < 0) {
13438865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
134447c6ae99SBarry Smith         }
134547c6ae99SBarry Smith         if (n26 >= 0) { /* right above */
1346ce00eea3SSatish Balay           x_t = lx[n26 % m];
134747c6ae99SBarry Smith           y_t = ly[(n26 % (m*n))/m];
134847c6ae99SBarry Smith           /* z_t = lz[n26 / (m*n)]; */
134947c6ae99SBarry Smith           s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
13508865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = s_t++;
1351ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ye-Ye < 0 && ze-Ze < 0) {
13528865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
135347c6ae99SBarry Smith         }
135447c6ae99SBarry Smith       }
135547c6ae99SBarry Smith     }
135647c6ae99SBarry Smith   }
135747c6ae99SBarry Smith   /*
135847c6ae99SBarry Smith      Set the local to global ordering in the global vector, this allows use
135947c6ae99SBarry Smith      of VecSetValuesLocal().
136047c6ae99SBarry Smith   */
1361*45b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr);
13623bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr);
136347c6ae99SBarry Smith 
1364ce00eea3SSatish Balay   ierr  = PetscFree2(bases,ldims);CHKERRQ(ierr);
1365ce00eea3SSatish Balay   dd->m = m;  dd->n  = n;  dd->p  = p;
1366ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
1367ce00eea3SSatish Balay   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = zs; dd->ze = ze;
1368ce00eea3SSatish Balay   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = Zs; dd->Ze = Ze;
1369ce00eea3SSatish Balay 
1370fcfd50ebSBarry Smith   ierr = VecDestroy(&local);CHKERRQ(ierr);
1371fcfd50ebSBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
1372ce00eea3SSatish Balay 
1373ce00eea3SSatish Balay   dd->gtol      = gtol;
1374ce00eea3SSatish Balay   dd->ltog      = ltog;
1375ce00eea3SSatish Balay   dd->base      = base;
1376ce00eea3SSatish Balay   da->ops->view = DMView_DA_3d;
13770298fd71SBarry Smith   dd->ltol      = NULL;
13780298fd71SBarry Smith   dd->ao        = NULL;
137947c6ae99SBarry Smith   PetscFunctionReturn(0);
138047c6ae99SBarry Smith }
1381564755cdSBarry Smith 
138247c6ae99SBarry Smith 
138347c6ae99SBarry Smith #undef __FUNCT__
1384aa219208SBarry Smith #define __FUNCT__ "DMDACreate3d"
138547c6ae99SBarry Smith /*@C
1386aa219208SBarry Smith    DMDACreate3d - Creates an object that will manage the communication of three-dimensional
138747c6ae99SBarry Smith    regular array data that is distributed across some processors.
138847c6ae99SBarry Smith 
138947c6ae99SBarry Smith    Collective on MPI_Comm
139047c6ae99SBarry Smith 
139147c6ae99SBarry Smith    Input Parameters:
139247c6ae99SBarry Smith +  comm - MPI communicator
13931321219cSEthan Coon .  bx,by,bz - type of ghost nodes the array have.
1394bff4a2f0SMatthew G. Knepley          Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC.
1395aa219208SBarry Smith .  stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX)
139647c6ae99SBarry 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
139747c6ae99SBarry Smith             from the command line with -da_grid_x <M> -da_grid_y <N> -da_grid_z <P>)
139847c6ae99SBarry Smith .  m,n,p - corresponding number of processors in each dimension
139947c6ae99SBarry Smith            (or PETSC_DECIDE to have calculated)
140047c6ae99SBarry Smith .  dof - number of degrees of freedom per node
140110d7c030SMatthew G Knepley .  s - stencil width
140210d7c030SMatthew G Knepley -  lx, ly, lz - arrays containing the number of nodes in each cell along
14030298fd71SBarry Smith           the x, y, and z coordinates, or NULL. If non-null, these
140447c6ae99SBarry Smith           must be of length as m,n,p and the corresponding
140547c6ae99SBarry Smith           m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of
140647c6ae99SBarry Smith           the ly[] must N, sum of the lz[] must be P
140747c6ae99SBarry Smith 
140847c6ae99SBarry Smith    Output Parameter:
140947c6ae99SBarry Smith .  da - the resulting distributed array object
141047c6ae99SBarry Smith 
141147c6ae99SBarry Smith    Options Database Key:
1412706a11cbSBarry Smith +  -dm_view - Calls DMView() at the conclusion of DMDACreate3d()
141347c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
141447c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
141547c6ae99SBarry Smith .  -da_grid_z <nz> - number of grid points in z direction, if P < 0
141647c6ae99SBarry Smith .  -da_processors_x <MX> - number of processors in x direction
141747c6ae99SBarry Smith .  -da_processors_y <MY> - number of processors in y direction
141847c6ae99SBarry Smith .  -da_processors_z <MZ> - number of processors in z direction
1419e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement ratio in x direction
1420e0f5d30fSBarry Smith .  -da_refine_y <ry> - refinement ratio in y direction
1421e0f5d30fSBarry Smith .  -da_refine_z <rz>- refinement ratio in z directio
1422e0f5d30fSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating it, , if M, N, or P < 0
142347c6ae99SBarry Smith 
142447c6ae99SBarry Smith    Level: beginner
142547c6ae99SBarry Smith 
142647c6ae99SBarry Smith    Notes:
1427aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
1428aa219208SBarry Smith    standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
142947c6ae99SBarry Smith    the standard 27-pt stencil.
143047c6ae99SBarry Smith 
1431aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
1432564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
1433564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
143447c6ae99SBarry Smith 
143547c6ae99SBarry Smith .keywords: distributed array, create, three-dimensional
143647c6ae99SBarry Smith 
1437aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate2d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
143899f0b487SRichard Tran Mills           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(),
1439d461ba97SBarry Smith           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
144047c6ae99SBarry Smith 
144147c6ae99SBarry Smith @*/
1442bff4a2f0SMatthew G. Knepley PetscErrorCode  DMDACreate3d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMBoundaryType bz,DMDAStencilType stencil_type,PetscInt M,
14439a42bb27SBarry 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)
144447c6ae99SBarry Smith {
144547c6ae99SBarry Smith   PetscErrorCode ierr;
144647c6ae99SBarry Smith 
144747c6ae99SBarry Smith   PetscFunctionBegin;
1448aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
1449aa219208SBarry Smith   ierr = DMDASetDim(*da, 3);CHKERRQ(ierr);
1450aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, P);CHKERRQ(ierr);
1451aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, p);CHKERRQ(ierr);
14521321219cSEthan Coon   ierr = DMDASetBoundaryType(*da, bx, by, bz);CHKERRQ(ierr);
1453aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
1454aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
1455aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
1456aa219208SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, lz);CHKERRQ(ierr);
145747c6ae99SBarry Smith   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
14589a42bb27SBarry Smith   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
14599a42bb27SBarry Smith   ierr = DMSetUp(*da);CHKERRQ(ierr);
146047c6ae99SBarry Smith   PetscFunctionReturn(0);
146147c6ae99SBarry Smith }
1462