xref: /petsc/src/dm/impls/da/da3.c (revision 3855c12bbba560cbb7c6a1785a9cd00b53ccb5e9)
17d0a6c19SBarry Smith 
247c6ae99SBarry Smith /*
347c6ae99SBarry Smith    Code for manipulating distributed regular 3d arrays in parallel.
447c6ae99SBarry Smith    File created by Peter Mell  7/14/95
547c6ae99SBarry Smith  */
647c6ae99SBarry Smith 
7c6db04a5SJed Brown #include <private/daimpl.h>     /*I   "petscdmda.h"    I*/
847c6ae99SBarry Smith 
947c6ae99SBarry Smith #undef __FUNCT__
109a42bb27SBarry Smith #define __FUNCT__ "DMView_DA_3d"
119a42bb27SBarry Smith PetscErrorCode DMView_DA_3d(DM da,PetscViewer viewer)
1247c6ae99SBarry Smith {
1347c6ae99SBarry Smith   PetscErrorCode ierr;
1447c6ae99SBarry Smith   PetscMPIInt    rank;
159a42bb27SBarry Smith   PetscBool      iascii,isdraw,isbinary;
1647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
179a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
189a42bb27SBarry Smith   PetscBool      ismatlab;
199a42bb27SBarry Smith #endif
2047c6ae99SBarry Smith 
2147c6ae99SBarry Smith   PetscFunctionBegin;
2247c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)da)->comm,&rank);CHKERRQ(ierr);
2347c6ae99SBarry Smith 
2447c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
2547c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
269a42bb27SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
279a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
289a42bb27SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr);
299a42bb27SBarry Smith #endif
3047c6ae99SBarry Smith   if (iascii) {
3147c6ae99SBarry Smith     PetscViewerFormat format;
3247c6ae99SBarry Smith 
337b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
3447c6ae99SBarry Smith     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
3547c6ae99SBarry Smith     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) {
36aa219208SBarry Smith       DMDALocalInfo info;
37aa219208SBarry Smith       ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
3847c6ae99SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %D N %D P %D m %D n %D p %D w %D s %D\n",rank,dd->M,dd->N,dd->P,dd->m,dd->n,dd->p,dd->w,dd->s);CHKERRQ(ierr);
3947c6ae99SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D, Y range of indices: %D %D, Z range of indices: %D %D\n",
4047c6ae99SBarry Smith                                                 info.xs,info.xs+info.xm,info.ys,info.ys+info.ym,info.zs,info.zs+info.zm);CHKERRQ(ierr);
4147c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX)
4247c6ae99SBarry Smith       if (dd->coordinates) {
4347c6ae99SBarry Smith         PetscInt        last;
4447c6ae99SBarry Smith         const PetscReal *coors;
4547c6ae99SBarry Smith         ierr = VecGetArrayRead(dd->coordinates,&coors);CHKERRQ(ierr);
4647c6ae99SBarry Smith         ierr = VecGetLocalSize(dd->coordinates,&last);CHKERRQ(ierr);
4747c6ae99SBarry Smith         last = last - 3;
4847c6ae99SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Lower left corner %G %G %G : Upper right %G %G %G\n",coors[0],coors[1],coors[2],coors[last],coors[last+1],coors[last+2]);CHKERRQ(ierr);
4947c6ae99SBarry Smith         ierr = VecRestoreArrayRead(dd->coordinates,&coors);CHKERRQ(ierr);
5047c6ae99SBarry Smith       }
5147c6ae99SBarry Smith #endif
5247c6ae99SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
537b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
54616157d6SJed Brown     } else {
55616157d6SJed Brown       ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr);
5647c6ae99SBarry Smith     }
5747c6ae99SBarry Smith   } else if (isdraw) {
5847c6ae99SBarry Smith     PetscDraw       draw;
5947c6ae99SBarry Smith     PetscReal     ymin = -1.0,ymax = (PetscReal)dd->N;
6047c6ae99SBarry Smith     PetscReal     xmin = -1.0,xmax = (PetscReal)((dd->M+2)*dd->P),x,y,ycoord,xcoord;
6147c6ae99SBarry Smith     PetscInt        k,plane,base,*idx;
6247c6ae99SBarry Smith     char       node[10];
6347c6ae99SBarry Smith     PetscBool  isnull;
6447c6ae99SBarry Smith 
6547c6ae99SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
6647c6ae99SBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
6747c6ae99SBarry Smith     ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
6847c6ae99SBarry Smith     ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
6947c6ae99SBarry Smith 
7047c6ae99SBarry Smith     /* first processor draw all node lines */
7147c6ae99SBarry Smith     if (!rank) {
7247c6ae99SBarry Smith       for (k=0; k<dd->P; k++) {
7347c6ae99SBarry Smith         ymin = 0.0; ymax = (PetscReal)(dd->N - 1);
7447c6ae99SBarry Smith         for (xmin=(PetscReal)(k*(dd->M+1)); xmin<(PetscReal)(dd->M+(k*(dd->M+1))); xmin++) {
7547c6ae99SBarry Smith           ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
7647c6ae99SBarry Smith         }
7747c6ae99SBarry Smith 
7847c6ae99SBarry Smith         xmin = (PetscReal)(k*(dd->M+1)); xmax = xmin + (PetscReal)(dd->M - 1);
7947c6ae99SBarry Smith         for (ymin=0; ymin<(PetscReal)dd->N; ymin++) {
8047c6ae99SBarry Smith           ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
8147c6ae99SBarry Smith         }
8247c6ae99SBarry Smith       }
8347c6ae99SBarry Smith     }
8447c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
8547c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
8647c6ae99SBarry Smith 
8747c6ae99SBarry Smith     for (k=0; k<dd->P; k++) {  /*Go through and draw for each plane*/
8847c6ae99SBarry Smith       if ((k >= dd->zs) && (k < dd->ze)) {
8947c6ae99SBarry Smith         /* draw my box */
9047c6ae99SBarry Smith         ymin = dd->ys;
9147c6ae99SBarry Smith         ymax = dd->ye - 1;
9247c6ae99SBarry Smith         xmin = dd->xs/dd->w    + (dd->M+1)*k;
9347c6ae99SBarry Smith         xmax =(dd->xe-1)/dd->w + (dd->M+1)*k;
9447c6ae99SBarry Smith 
9547c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
9647c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
9747c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
9847c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
9947c6ae99SBarry Smith 
10047c6ae99SBarry Smith         xmin = dd->xs/dd->w;
10147c6ae99SBarry Smith         xmax =(dd->xe-1)/dd->w;
10247c6ae99SBarry Smith 
10347c6ae99SBarry Smith         /* put in numbers*/
10447c6ae99SBarry Smith         base = (dd->base+(dd->xe-dd->xs)*(dd->ye-dd->ys)*(k-dd->zs))/dd->w;
10547c6ae99SBarry Smith 
10647c6ae99SBarry Smith         /* Identify which processor owns the box */
10747c6ae99SBarry Smith         sprintf(node,"%d",rank);
10847c6ae99SBarry Smith         ierr = PetscDrawString(draw,xmin+(dd->M+1)*k+.2,ymin+.3,PETSC_DRAW_RED,node);CHKERRQ(ierr);
10947c6ae99SBarry Smith 
11047c6ae99SBarry Smith         for (y=ymin; y<=ymax; y++) {
11147c6ae99SBarry Smith           for (x=xmin+(dd->M+1)*k; x<=xmax+(dd->M+1)*k; x++) {
11247c6ae99SBarry Smith             sprintf(node,"%d",(int)base++);
11347c6ae99SBarry Smith             ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
11447c6ae99SBarry Smith           }
11547c6ae99SBarry Smith         }
11647c6ae99SBarry Smith 
11747c6ae99SBarry Smith       }
11847c6ae99SBarry Smith     }
11947c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
12047c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
12147c6ae99SBarry Smith 
12247c6ae99SBarry Smith     for (k=0-dd->s; k<dd->P+dd->s; k++) {
12347c6ae99SBarry Smith       /* Go through and draw for each plane */
12447c6ae99SBarry Smith       if ((k >= dd->Zs) && (k < dd->Ze)) {
12547c6ae99SBarry Smith 
12647c6ae99SBarry Smith         /* overlay ghost numbers, useful for error checking */
12747c6ae99SBarry Smith         base = (dd->Xe-dd->Xs)*(dd->Ye-dd->Ys)*(k-dd->Zs); idx = dd->idx;
12847c6ae99SBarry Smith         plane=k;
12947c6ae99SBarry Smith         /* Keep z wrap around points on the dradrawg */
13047c6ae99SBarry Smith         if (k<0)    { plane=dd->P+k; }
13147c6ae99SBarry Smith         if (k>=dd->P) { plane=k-dd->P; }
13247c6ae99SBarry Smith         ymin = dd->Ys; ymax = dd->Ye;
13347c6ae99SBarry Smith         xmin = (dd->M+1)*plane*dd->w;
13447c6ae99SBarry Smith         xmax = (dd->M+1)*plane*dd->w+dd->M*dd->w;
13547c6ae99SBarry Smith         for (y=ymin; y<ymax; y++) {
13647c6ae99SBarry Smith           for (x=xmin+dd->Xs; x<xmin+dd->Xe; x+=dd->w) {
13747c6ae99SBarry Smith             sprintf(node,"%d",(int)(idx[base]/dd->w));
13847c6ae99SBarry Smith             ycoord = y;
13947c6ae99SBarry Smith             /*Keep y wrap around points on drawing */
14047c6ae99SBarry Smith             if (y<0)      { ycoord = dd->N+y; }
14147c6ae99SBarry Smith 
14247c6ae99SBarry Smith             if (y>=dd->N) { ycoord = y-dd->N; }
14347c6ae99SBarry Smith             xcoord = x;   /* Keep x wrap points on drawing */
14447c6ae99SBarry Smith 
14547c6ae99SBarry Smith             if (x<xmin)  { xcoord = xmax - (xmin-x); }
14647c6ae99SBarry Smith             if (x>=xmax) { xcoord = xmin + (x-xmax); }
14747c6ae99SBarry Smith             ierr = PetscDrawString(draw,xcoord/dd->w,ycoord,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
14847c6ae99SBarry Smith             base+=dd->w;
14947c6ae99SBarry Smith           }
15047c6ae99SBarry Smith         }
15147c6ae99SBarry Smith       }
15247c6ae99SBarry Smith     }
15347c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
15447c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
1559a42bb27SBarry Smith   } else if (isbinary){
1569a42bb27SBarry Smith     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
1579a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1589a42bb27SBarry Smith   } else if (ismatlab) {
1599a42bb27SBarry Smith     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
1609a42bb27SBarry Smith #endif
161aa219208SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for DMDA 1d",((PetscObject)viewer)->type_name);
16247c6ae99SBarry Smith   PetscFunctionReturn(0);
16347c6ae99SBarry Smith }
16447c6ae99SBarry Smith 
16547c6ae99SBarry Smith #undef __FUNCT__
1669a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_3D"
1677087cfbeSBarry Smith PetscErrorCode  DMSetUp_DA_3D(DM da)
16847c6ae99SBarry Smith {
16947c6ae99SBarry Smith   DM_DA                  *dd           = (DM_DA*)da->data;
17047c6ae99SBarry Smith   const PetscInt         M            = dd->M;
17147c6ae99SBarry Smith   const PetscInt         N            = dd->N;
17247c6ae99SBarry Smith   const PetscInt         P            = dd->P;
17347c6ae99SBarry Smith   PetscInt               m            = dd->m;
17447c6ae99SBarry Smith   PetscInt               n            = dd->n;
17547c6ae99SBarry Smith   PetscInt               p            = dd->p;
17647c6ae99SBarry Smith   const PetscInt         dof          = dd->w;
17747c6ae99SBarry Smith   const PetscInt         s            = dd->s;
1781321219cSEthan Coon   const DMDABoundaryType bx         = dd->bx;
1791321219cSEthan Coon   const DMDABoundaryType by         = dd->by;
1801321219cSEthan Coon   const DMDABoundaryType bz         = dd->bz;
181aa219208SBarry Smith   const DMDAStencilType  stencil_type = dd->stencil_type;
18247c6ae99SBarry Smith   PetscInt               *lx           = dd->lx;
18347c6ae99SBarry Smith   PetscInt               *ly           = dd->ly;
18447c6ae99SBarry Smith   PetscInt               *lz           = dd->lz;
18547c6ae99SBarry Smith   MPI_Comm               comm;
18647c6ae99SBarry Smith   PetscMPIInt            rank,size;
187ce00eea3SSatish Balay   PetscInt               xs = 0,xe,ys = 0,ye,zs = 0,ze,x = 0,y = 0,z = 0;
188ce00eea3SSatish Balay   PetscInt               Xs,Xe,Ys,Ye,Zs,Ze,IXs,IXe,IYs,IYe,IZs,IZe,start,end,pm;
189ce00eea3SSatish Balay   PetscInt               left,right,up,down,bottom,top,i,j,k,*idx,*idx_cpy,nn;
190ce00eea3SSatish Balay   const PetscInt         *idx_full;
19147c6ae99SBarry Smith   PetscInt               n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14;
19247c6ae99SBarry Smith   PetscInt               n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26;
193db87c5ecSEthan Coon   PetscInt               *bases,*ldims,base,x_t,y_t,z_t,s_t,count,s_x,s_y,s_z;
19447c6ae99SBarry Smith   PetscInt               sn0 = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0;
19547c6ae99SBarry Smith   PetscInt               sn8 = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0;
19647c6ae99SBarry Smith   PetscInt               sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0;
19747c6ae99SBarry Smith   Vec                    local,global;
19847c6ae99SBarry Smith   VecScatter             ltog,gtol;
199ce00eea3SSatish Balay   IS                     to,from,ltogis;
20047c6ae99SBarry Smith   PetscErrorCode         ierr;
20147c6ae99SBarry Smith 
20247c6ae99SBarry Smith   PetscFunctionBegin;
20347c6ae99SBarry Smith   if (dof < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
20447c6ae99SBarry Smith   if (s < 0) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
20547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr);
206*3855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES)
207*3855c12bSBarry 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);
208*3855c12bSBarry Smith #endif
209*3855c12bSBarry Smith 
21047c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
21147c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
21247c6ae99SBarry Smith 
21347c6ae99SBarry Smith   dd->dim = 3;
21447c6ae99SBarry Smith   ierr = PetscMalloc(dof*sizeof(char*),&dd->fieldname);CHKERRQ(ierr);
21547c6ae99SBarry Smith   ierr = PetscMemzero(dd->fieldname,dof*sizeof(char*));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 */
2408f1a2a5eSBarry Smith     m = (int)(0.5 + sqrt(((double)M)*((double)size)/((double)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 */
2518f1a2a5eSBarry Smith     m = (int)(0.5 + sqrt(((double)M)*((double)size)/((double)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 */
2628f1a2a5eSBarry Smith     n = (int)(0.5 + sqrt(((double)N)*((double)size)/((double)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 */
2738f1a2a5eSBarry Smith     n = (PetscInt)(0.5 + pow(((double)N*N)*((double)size)/((double)P*M),(double)(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;
2818f1a2a5eSBarry Smith     m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)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;}
28947c6ae99SBarry Smith   } else if (m*n*p != size) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
29047c6ae99SBarry Smith 
29147c6ae99SBarry Smith   if (m*n*p != size) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_PLIB,"Could not find good partition");
29247c6ae99SBarry Smith   if (M < m) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
29347c6ae99SBarry Smith   if (N < n) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
29447c6ae99SBarry Smith   if (P < p) SETERRQ2(((PetscObject)da)->comm,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) {
30247c6ae99SBarry Smith     ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
30347c6ae99SBarry Smith     lx = dd->lx;
30447c6ae99SBarry Smith     for (i=0; i<m; i++) {
30547c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > (i % m));
30647c6ae99SBarry Smith     }
30747c6ae99SBarry Smith   }
30847c6ae99SBarry Smith   x  = lx[rank % m];
30947c6ae99SBarry Smith   xs = 0;
31047c6ae99SBarry Smith   for (i=0; i<(rank%m); i++) { xs += lx[i];}
311bcea557cSEthan Coon   if ((x < s) && ((m > 1) || (bx == DMDA_BOUNDARY_PERIODIC))) {
312bcea557cSEthan Coon     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s);
313bcea557cSEthan Coon   }
31447c6ae99SBarry Smith 
31547c6ae99SBarry Smith   if (!ly) {
31647c6ae99SBarry Smith     ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
31747c6ae99SBarry Smith     ly = dd->ly;
31847c6ae99SBarry Smith     for (i=0; i<n; i++) {
31947c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > (i % n));
32047c6ae99SBarry Smith     }
32147c6ae99SBarry Smith   }
32247c6ae99SBarry Smith   y  = ly[(rank % (m*n))/m];
323bcea557cSEthan Coon   if ((y < s) && ((n > 1) || (by == DMDA_BOUNDARY_PERIODIC))) {
324bcea557cSEthan Coon     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s);
325bcea557cSEthan Coon   }
32647c6ae99SBarry Smith   ys = 0;
32747c6ae99SBarry Smith   for (i=0; i<(rank % (m*n))/m; i++) { ys += ly[i];}
32847c6ae99SBarry Smith 
32947c6ae99SBarry Smith   if (!lz) {
33047c6ae99SBarry Smith     ierr = PetscMalloc(p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr);
33147c6ae99SBarry Smith     lz = dd->lz;
33247c6ae99SBarry Smith     for (i=0; i<p; i++) {
33347c6ae99SBarry Smith       lz[i] = P/p + ((P % p) > (i % p));
33447c6ae99SBarry Smith     }
33547c6ae99SBarry Smith   }
33647c6ae99SBarry Smith   z  = lz[rank/(m*n)];
337bcea557cSEthan Coon 
338fdc81ce8SEthan Coon   /* note this is different than x- and y-, as we will handle as an important special
339fdc81ce8SEthan Coon    case when p=P=1 and DMDA_BOUNDARY_PERIODIC and s > z.  This is to deal with 2D problems
340fdc81ce8SEthan Coon    in a 3D code.  Additional code for this case is noted with "2d case" comments */
341fdc81ce8SEthan Coon   if ((z < s) && ((p > 1) || ((P > 1) && bz == DMDA_BOUNDARY_PERIODIC))) {
342bcea557cSEthan Coon     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local z-width of domain z %D is smaller than stencil width s %D",z,s);
343bcea557cSEthan Coon   }
34447c6ae99SBarry Smith   zs = 0;
34547c6ae99SBarry Smith   for (i=0; i<(rank/(m*n)); i++) { zs += lz[i];}
34647c6ae99SBarry Smith   ye = ys + y;
34747c6ae99SBarry Smith   xe = xs + x;
34847c6ae99SBarry Smith   ze = zs + z;
34947c6ae99SBarry Smith 
35047c6ae99SBarry Smith   /* determine ghost region */
35147c6ae99SBarry Smith   /* Assume No Periodicity */
352ce00eea3SSatish Balay   if (xs-s > 0) { Xs = xs - s; IXs = xs - s; } else { Xs = 0; IXs = 0; }
353ce00eea3SSatish Balay   if (xe+s <= M) { Xe = xe + s; IXe = xe + s; } else { Xe = M; IXe = M; }
354ce00eea3SSatish Balay   if (ys-s > 0) { Ys = ys - s; IYs = ys - s; } else { Ys = 0; IYs = 0; }
355ce00eea3SSatish Balay   if (ye+s <= N) { Ye = ye + s; IYe = ye + s; } else { Ye = N; IYe = N; }
356ce00eea3SSatish Balay   if (zs-s > 0) { Zs = zs - s; IZs = zs - s; } else { Zs = 0; IZs = 0; }
357ce00eea3SSatish Balay   if (ze+s <= P) { Ze = ze + s; IZe = ze + s; } else { Ze = P; IZe = P; }
35847c6ae99SBarry Smith 
359ce00eea3SSatish Balay   /* fix for periodicity/ghosted */
3601321219cSEthan Coon   if (bx) { Xs = xs - s; Xe = xe + s; }
3611321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC) { IXs = xs - s; IXe = xe + s; }
3621321219cSEthan Coon   if (by) { Ys = ys - s; Ye = ye + s; }
3631321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC) { IYs = ys - s; IYe = ye + s; }
3641321219cSEthan Coon   if (bz) { Zs = zs - s; Ze = ze + s; }
3651321219cSEthan Coon   if (bz == DMDA_BOUNDARY_PERIODIC) { IZs = zs - s; IZe = ze + s; }
36647c6ae99SBarry Smith 
36747c6ae99SBarry Smith   /* Resize all X parameters to reflect w */
368ce00eea3SSatish Balay   s_x = s;
36947c6ae99SBarry Smith   s_y  = s;
37047c6ae99SBarry Smith   s_z  = s;
37147c6ae99SBarry Smith 
37247c6ae99SBarry Smith   /* determine starting point of each processor */
37347c6ae99SBarry Smith   nn       = x*y*z;
37447c6ae99SBarry Smith   ierr     = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr);
37547c6ae99SBarry Smith   ierr     = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
37647c6ae99SBarry Smith   bases[0] = 0;
37747c6ae99SBarry Smith   for (i=1; i<=size; i++) {
37847c6ae99SBarry Smith     bases[i] = ldims[i-1];
37947c6ae99SBarry Smith   }
38047c6ae99SBarry Smith   for (i=1; i<=size; i++) {
38147c6ae99SBarry Smith     bases[i] += bases[i-1];
38247c6ae99SBarry Smith   }
383ce00eea3SSatish Balay   base = bases[rank]*dof;
38447c6ae99SBarry Smith 
38547c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
386ce00eea3SSatish Balay   dd->Nlocal = x*y*z*dof;
38747c6ae99SBarry Smith   ierr = VecCreateMPIWithArray(comm,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr);
38847c6ae99SBarry Smith   ierr = VecSetBlockSize(global,dof);CHKERRQ(ierr);
389ce00eea3SSatish Balay   dd->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs)*dof;
390ce00eea3SSatish Balay   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dd->nlocal,0,&local);CHKERRQ(ierr);
39147c6ae99SBarry Smith   ierr = VecSetBlockSize(local,dof);CHKERRQ(ierr);
39247c6ae99SBarry Smith 
39347c6ae99SBarry Smith   /* generate appropriate vector scatters */
39447c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
39547c6ae99SBarry Smith   ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr);
396ce00eea3SSatish Balay   ierr = ISCreateStride(comm,x*y*z*dof,start,1,&to);CHKERRQ(ierr);
39747c6ae99SBarry Smith 
398db87c5ecSEthan Coon   count = x*y*z;
399ce00eea3SSatish Balay   ierr = PetscMalloc(x*y*z*sizeof(PetscInt),&idx);CHKERRQ(ierr);
400ce00eea3SSatish Balay   left   = xs - Xs; right = left + x;
40147c6ae99SBarry Smith   bottom = ys - Ys; top = bottom + y;
40247c6ae99SBarry Smith   down   = zs - Zs; up  = down + z;
40347c6ae99SBarry Smith   count  = 0;
40447c6ae99SBarry Smith   for (i=down; i<up; i++) {
40547c6ae99SBarry Smith     for (j=bottom; j<top; j++) {
406ce00eea3SSatish Balay       for (k=left; k<right; k++) {
407ce00eea3SSatish Balay         idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
40847c6ae99SBarry Smith       }
40947c6ae99SBarry Smith     }
41047c6ae99SBarry Smith   }
41147c6ae99SBarry Smith 
41247c6ae99SBarry Smith   ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
41347c6ae99SBarry Smith   ierr = VecScatterCreate(local,from,global,to,&ltog);CHKERRQ(ierr);
41447c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,ltog);CHKERRQ(ierr);
415fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
416fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
41747c6ae99SBarry Smith 
418ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
419ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
420aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_BOX) {
421db87c5ecSEthan Coon     count = (IXe-IXs)*(IYe-IYs)*(IZe-IZs);
422db87c5ecSEthan Coon     ierr  = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
423ce00eea3SSatish Balay 
424ce00eea3SSatish Balay     left   = IXs - Xs; right = left + (IXe-IXs);
425ce00eea3SSatish Balay     bottom = IYs - Ys; top = bottom + (IYe-IYs);
426ce00eea3SSatish Balay     down   = IZs - Zs; up  = down + (IZe-IZs);
427ce00eea3SSatish Balay     count = 0;
428ce00eea3SSatish Balay     for (i=down; i<up; i++) {
429ce00eea3SSatish Balay       for (j=bottom; j<top; j++) {
430ce00eea3SSatish Balay         for (k=left; k<right; k++) {
431ce00eea3SSatish Balay           idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
432ce00eea3SSatish Balay         }
433ce00eea3SSatish Balay       }
434ce00eea3SSatish Balay     }
435ce00eea3SSatish Balay     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
436ce00eea3SSatish Balay 
43747c6ae99SBarry Smith   } else {
43847c6ae99SBarry Smith     /* This is way ugly! We need to list the funny cross type region */
439db87c5ecSEthan Coon     count = ((ys-IYs) + (IYe-ye))*x*z + ((xs-IXs) + (IXe-xe))*y*z + ((zs-IZs) + (IZe-ze))*x*y + x*y*z;
440db87c5ecSEthan Coon     ierr   = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
441ce00eea3SSatish Balay 
442ce00eea3SSatish Balay     left   = xs - Xs; right = left + x;
44347c6ae99SBarry Smith     bottom = ys - Ys; top = bottom + y;
44447c6ae99SBarry Smith     down   = zs - Zs;   up  = down + z;
44547c6ae99SBarry Smith     count  = 0;
446ce00eea3SSatish Balay     /* the bottom chunck */
447ce00eea3SSatish Balay     for (i=(IZs-Zs); i<down; i++) {
44847c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
449ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
45047c6ae99SBarry Smith       }
45147c6ae99SBarry Smith     }
45247c6ae99SBarry Smith     /* the middle piece */
45347c6ae99SBarry Smith     for (i=down; i<up; i++) {
45447c6ae99SBarry Smith       /* front */
455ce00eea3SSatish Balay       for (j=(IYs-Ys); j<bottom; j++) {
456ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
45747c6ae99SBarry Smith       }
45847c6ae99SBarry Smith       /* middle */
45947c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
460ce00eea3SSatish Balay         for (k=IXs-Xs; k<IXe-Xs; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
46147c6ae99SBarry Smith       }
46247c6ae99SBarry Smith       /* back */
463ce00eea3SSatish Balay       for (j=top; j<top+IYe-ye; j++) {
464ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
46547c6ae99SBarry Smith       }
46647c6ae99SBarry Smith     }
46747c6ae99SBarry Smith     /* the top piece */
468ce00eea3SSatish Balay     for (i=up; i<up+IZe-ze; i++) {
46947c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
470ce00eea3SSatish Balay         for (k=left; k<right; k++) idx[count++] = (i*(Ye-Ys) + j)*(Xe-Xs) + k;
47147c6ae99SBarry Smith       }
47247c6ae99SBarry Smith     }
47347c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
47447c6ae99SBarry Smith   }
47547c6ae99SBarry Smith 
47647c6ae99SBarry Smith   /* determine who lies on each side of use stored in    n24 n25 n26
47747c6ae99SBarry Smith                                                          n21 n22 n23
47847c6ae99SBarry Smith                                                          n18 n19 n20
47947c6ae99SBarry Smith 
48047c6ae99SBarry Smith                                                          n15 n16 n17
48147c6ae99SBarry Smith                                                          n12     n14
48247c6ae99SBarry Smith                                                          n9  n10 n11
48347c6ae99SBarry Smith 
48447c6ae99SBarry Smith                                                          n6  n7  n8
48547c6ae99SBarry Smith                                                          n3  n4  n5
48647c6ae99SBarry Smith                                                          n0  n1  n2
48747c6ae99SBarry Smith   */
48847c6ae99SBarry Smith 
48947c6ae99SBarry Smith   /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
49047c6ae99SBarry Smith   /* Assume Nodes are Internal to the Cube */
49147c6ae99SBarry Smith   n0  = rank - m*n - m - 1;
49247c6ae99SBarry Smith   n1  = rank - m*n - m;
49347c6ae99SBarry Smith   n2  = rank - m*n - m + 1;
49447c6ae99SBarry Smith   n3  = rank - m*n -1;
49547c6ae99SBarry Smith   n4  = rank - m*n;
49647c6ae99SBarry Smith   n5  = rank - m*n + 1;
49747c6ae99SBarry Smith   n6  = rank - m*n + m - 1;
49847c6ae99SBarry Smith   n7  = rank - m*n + m;
49947c6ae99SBarry Smith   n8  = rank - m*n + m + 1;
50047c6ae99SBarry Smith 
50147c6ae99SBarry Smith   n9  = rank - m - 1;
50247c6ae99SBarry Smith   n10 = rank - m;
50347c6ae99SBarry Smith   n11 = rank - m + 1;
50447c6ae99SBarry Smith   n12 = rank - 1;
50547c6ae99SBarry Smith   n14 = rank + 1;
50647c6ae99SBarry Smith   n15 = rank + m - 1;
50747c6ae99SBarry Smith   n16 = rank + m;
50847c6ae99SBarry Smith   n17 = rank + m + 1;
50947c6ae99SBarry Smith 
51047c6ae99SBarry Smith   n18 = rank + m*n - m - 1;
51147c6ae99SBarry Smith   n19 = rank + m*n - m;
51247c6ae99SBarry Smith   n20 = rank + m*n - m + 1;
51347c6ae99SBarry Smith   n21 = rank + m*n - 1;
51447c6ae99SBarry Smith   n22 = rank + m*n;
51547c6ae99SBarry Smith   n23 = rank + m*n + 1;
51647c6ae99SBarry Smith   n24 = rank + m*n + m - 1;
51747c6ae99SBarry Smith   n25 = rank + m*n + m;
51847c6ae99SBarry Smith   n26 = rank + m*n + m + 1;
51947c6ae99SBarry Smith 
52047c6ae99SBarry Smith   /* Assume Pieces are on Faces of Cube */
52147c6ae99SBarry Smith 
52247c6ae99SBarry Smith   if (xs == 0) { /* First assume not corner or edge */
52347c6ae99SBarry Smith     n0  = rank       -1 - (m*n);
52447c6ae99SBarry Smith     n3  = rank + m   -1 - (m*n);
52547c6ae99SBarry Smith     n6  = rank + 2*m -1 - (m*n);
52647c6ae99SBarry Smith     n9  = rank       -1;
52747c6ae99SBarry Smith     n12 = rank + m   -1;
52847c6ae99SBarry Smith     n15 = rank + 2*m -1;
52947c6ae99SBarry Smith     n18 = rank       -1 + (m*n);
53047c6ae99SBarry Smith     n21 = rank + m   -1 + (m*n);
53147c6ae99SBarry Smith     n24 = rank + 2*m -1 + (m*n);
53247c6ae99SBarry Smith    }
53347c6ae99SBarry Smith 
534ce00eea3SSatish Balay   if (xe == M) { /* First assume not corner or edge */
53547c6ae99SBarry Smith     n2  = rank -2*m +1 - (m*n);
53647c6ae99SBarry Smith     n5  = rank - m  +1 - (m*n);
53747c6ae99SBarry Smith     n8  = rank      +1 - (m*n);
53847c6ae99SBarry Smith     n11 = rank -2*m +1;
53947c6ae99SBarry Smith     n14 = rank - m  +1;
54047c6ae99SBarry Smith     n17 = rank      +1;
54147c6ae99SBarry Smith     n20 = rank -2*m +1 + (m*n);
54247c6ae99SBarry Smith     n23 = rank - m  +1 + (m*n);
54347c6ae99SBarry Smith     n26 = rank      +1 + (m*n);
54447c6ae99SBarry Smith   }
54547c6ae99SBarry Smith 
54647c6ae99SBarry Smith   if (ys==0) { /* First assume not corner or edge */
54747c6ae99SBarry Smith     n0  = rank + m * (n-1) -1 - (m*n);
54847c6ae99SBarry Smith     n1  = rank + m * (n-1)    - (m*n);
54947c6ae99SBarry Smith     n2  = rank + m * (n-1) +1 - (m*n);
55047c6ae99SBarry Smith     n9  = rank + m * (n-1) -1;
55147c6ae99SBarry Smith     n10 = rank + m * (n-1);
55247c6ae99SBarry Smith     n11 = rank + m * (n-1) +1;
55347c6ae99SBarry Smith     n18 = rank + m * (n-1) -1 + (m*n);
55447c6ae99SBarry Smith     n19 = rank + m * (n-1)    + (m*n);
55547c6ae99SBarry Smith     n20 = rank + m * (n-1) +1 + (m*n);
55647c6ae99SBarry Smith   }
55747c6ae99SBarry Smith 
55847c6ae99SBarry Smith   if (ye == N) { /* First assume not corner or edge */
55947c6ae99SBarry Smith     n6  = rank - m * (n-1) -1 - (m*n);
56047c6ae99SBarry Smith     n7  = rank - m * (n-1)    - (m*n);
56147c6ae99SBarry Smith     n8  = rank - m * (n-1) +1 - (m*n);
56247c6ae99SBarry Smith     n15 = rank - m * (n-1) -1;
56347c6ae99SBarry Smith     n16 = rank - m * (n-1);
56447c6ae99SBarry Smith     n17 = rank - m * (n-1) +1;
56547c6ae99SBarry Smith     n24 = rank - m * (n-1) -1 + (m*n);
56647c6ae99SBarry Smith     n25 = rank - m * (n-1)    + (m*n);
56747c6ae99SBarry Smith     n26 = rank - m * (n-1) +1 + (m*n);
56847c6ae99SBarry Smith   }
56947c6ae99SBarry Smith 
57047c6ae99SBarry Smith   if (zs == 0) { /* First assume not corner or edge */
57147c6ae99SBarry Smith     n0 = size - (m*n) + rank - m - 1;
57247c6ae99SBarry Smith     n1 = size - (m*n) + rank - m;
57347c6ae99SBarry Smith     n2 = size - (m*n) + rank - m + 1;
57447c6ae99SBarry Smith     n3 = size - (m*n) + rank - 1;
57547c6ae99SBarry Smith     n4 = size - (m*n) + rank;
57647c6ae99SBarry Smith     n5 = size - (m*n) + rank + 1;
57747c6ae99SBarry Smith     n6 = size - (m*n) + rank + m - 1;
57847c6ae99SBarry Smith     n7 = size - (m*n) + rank + m ;
57947c6ae99SBarry Smith     n8 = size - (m*n) + rank + m + 1;
58047c6ae99SBarry Smith   }
58147c6ae99SBarry Smith 
58247c6ae99SBarry Smith   if (ze == P) { /* First assume not corner or edge */
58347c6ae99SBarry Smith     n18 = (m*n) - (size-rank) - m - 1;
58447c6ae99SBarry Smith     n19 = (m*n) - (size-rank) - m;
58547c6ae99SBarry Smith     n20 = (m*n) - (size-rank) - m + 1;
58647c6ae99SBarry Smith     n21 = (m*n) - (size-rank) - 1;
58747c6ae99SBarry Smith     n22 = (m*n) - (size-rank);
58847c6ae99SBarry Smith     n23 = (m*n) - (size-rank) + 1;
58947c6ae99SBarry Smith     n24 = (m*n) - (size-rank) + m - 1;
59047c6ae99SBarry Smith     n25 = (m*n) - (size-rank) + m;
59147c6ae99SBarry Smith     n26 = (m*n) - (size-rank) + m + 1;
59247c6ae99SBarry Smith   }
59347c6ae99SBarry Smith 
59447c6ae99SBarry Smith   if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */
59547c6ae99SBarry Smith     n0 = size - m*n + rank + m-1 - m;
59647c6ae99SBarry Smith     n3 = size - m*n + rank + m-1;
59747c6ae99SBarry Smith     n6 = size - m*n + rank + m-1 + m;
59847c6ae99SBarry Smith   }
59947c6ae99SBarry Smith 
60047c6ae99SBarry Smith   if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */
60147c6ae99SBarry Smith     n18 = m*n - (size - rank) + m-1 - m;
60247c6ae99SBarry Smith     n21 = m*n - (size - rank) + m-1;
60347c6ae99SBarry Smith     n24 = m*n - (size - rank) + m-1 + m;
60447c6ae99SBarry Smith   }
60547c6ae99SBarry Smith 
60647c6ae99SBarry Smith   if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */
60747c6ae99SBarry Smith     n0  = rank + m*n -1 - m*n;
60847c6ae99SBarry Smith     n9  = rank + m*n -1;
60947c6ae99SBarry Smith     n18 = rank + m*n -1 + m*n;
61047c6ae99SBarry Smith   }
61147c6ae99SBarry Smith 
61247c6ae99SBarry Smith   if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */
61347c6ae99SBarry Smith     n6  = rank - m*(n-1) + m-1 - m*n;
61447c6ae99SBarry Smith     n15 = rank - m*(n-1) + m-1;
61547c6ae99SBarry Smith     n24 = rank - m*(n-1) + m-1 + m*n;
61647c6ae99SBarry Smith   }
61747c6ae99SBarry Smith 
618ce00eea3SSatish Balay   if ((xe==M) && (zs==0)) { /* Assume an edge, not corner */
61947c6ae99SBarry Smith     n2 = size - (m*n-rank) - (m-1) - m;
62047c6ae99SBarry Smith     n5 = size - (m*n-rank) - (m-1);
62147c6ae99SBarry Smith     n8 = size - (m*n-rank) - (m-1) + m;
62247c6ae99SBarry Smith   }
62347c6ae99SBarry Smith 
624ce00eea3SSatish Balay   if ((xe==M) && (ze==P)) { /* Assume an edge, not corner */
62547c6ae99SBarry Smith     n20 = m*n - (size - rank) - (m-1) - m;
62647c6ae99SBarry Smith     n23 = m*n - (size - rank) - (m-1);
62747c6ae99SBarry Smith     n26 = m*n - (size - rank) - (m-1) + m;
62847c6ae99SBarry Smith   }
62947c6ae99SBarry Smith 
630ce00eea3SSatish Balay   if ((xe==M) && (ys==0)) { /* Assume an edge, not corner */
63147c6ae99SBarry Smith     n2  = rank + m*(n-1) - (m-1) - m*n;
63247c6ae99SBarry Smith     n11 = rank + m*(n-1) - (m-1);
63347c6ae99SBarry Smith     n20 = rank + m*(n-1) - (m-1) + m*n;
63447c6ae99SBarry Smith   }
63547c6ae99SBarry Smith 
636ce00eea3SSatish Balay   if ((xe==M) && (ye==N)) { /* Assume an edge, not corner */
63747c6ae99SBarry Smith     n8  = rank - m*n +1 - m*n;
63847c6ae99SBarry Smith     n17 = rank - m*n +1;
63947c6ae99SBarry Smith     n26 = rank - m*n +1 + m*n;
64047c6ae99SBarry Smith   }
64147c6ae99SBarry Smith 
64247c6ae99SBarry Smith   if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */
64347c6ae99SBarry Smith     n0 = size - m + rank -1;
64447c6ae99SBarry Smith     n1 = size - m + rank;
64547c6ae99SBarry Smith     n2 = size - m + rank +1;
64647c6ae99SBarry Smith   }
64747c6ae99SBarry Smith 
64847c6ae99SBarry Smith   if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */
64947c6ae99SBarry Smith     n18 = m*n - (size - rank) + m*(n-1) -1;
65047c6ae99SBarry Smith     n19 = m*n - (size - rank) + m*(n-1);
65147c6ae99SBarry Smith     n20 = m*n - (size - rank) + m*(n-1) +1;
65247c6ae99SBarry Smith   }
65347c6ae99SBarry Smith 
65447c6ae99SBarry Smith   if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */
65547c6ae99SBarry Smith     n6 = size - (m*n-rank) - m * (n-1) -1;
65647c6ae99SBarry Smith     n7 = size - (m*n-rank) - m * (n-1);
65747c6ae99SBarry Smith     n8 = size - (m*n-rank) - m * (n-1) +1;
65847c6ae99SBarry Smith   }
65947c6ae99SBarry Smith 
66047c6ae99SBarry Smith   if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */
66147c6ae99SBarry Smith     n24 = rank - (size-m) -1;
66247c6ae99SBarry Smith     n25 = rank - (size-m);
66347c6ae99SBarry Smith     n26 = rank - (size-m) +1;
66447c6ae99SBarry Smith   }
66547c6ae99SBarry Smith 
66647c6ae99SBarry Smith   /* Check for Corners */
66747c6ae99SBarry Smith   if ((xs==0)   && (ys==0) && (zs==0)) { n0  = size -1;}
66847c6ae99SBarry Smith   if ((xs==0)   && (ys==0) && (ze==P)) { n18 = m*n-1;}
66947c6ae99SBarry Smith   if ((xs==0)   && (ye==N) && (zs==0)) { n6  = (size-1)-m*(n-1);}
67047c6ae99SBarry Smith   if ((xs==0)   && (ye==N) && (ze==P)) { n24 = m-1;}
671ce00eea3SSatish Balay   if ((xe==M) && (ys==0) && (zs==0)) { n2  = size-m;}
672ce00eea3SSatish Balay   if ((xe==M) && (ys==0) && (ze==P)) { n20 = m*n-m;}
673ce00eea3SSatish Balay   if ((xe==M) && (ye==N) && (zs==0)) { n8  = size-m*n;}
674ce00eea3SSatish Balay   if ((xe==M) && (ye==N) && (ze==P)) { n26 = 0;}
67547c6ae99SBarry Smith 
67647c6ae99SBarry Smith   /* Check for when not X,Y, and Z Periodic */
67747c6ae99SBarry Smith 
67847c6ae99SBarry Smith   /* If not X periodic */
6791321219cSEthan Coon   if (bx != DMDA_BOUNDARY_PERIODIC) {
68047c6ae99SBarry Smith     if (xs==0)   {n0  = n3  = n6  = n9  = n12 = n15 = n18 = n21 = n24 = -2;}
681ce00eea3SSatish Balay     if (xe==M) {n2  = n5  = n8  = n11 = n14 = n17 = n20 = n23 = n26 = -2;}
68247c6ae99SBarry Smith   }
68347c6ae99SBarry Smith 
68447c6ae99SBarry Smith   /* If not Y periodic */
6851321219cSEthan Coon   if (by != DMDA_BOUNDARY_PERIODIC) {
68647c6ae99SBarry Smith     if (ys==0)   {n0  = n1  = n2  = n9  = n10 = n11 = n18 = n19 = n20 = -2;}
68747c6ae99SBarry Smith     if (ye==N)   {n6  = n7  = n8  = n15 = n16 = n17 = n24 = n25 = n26 = -2;}
68847c6ae99SBarry Smith   }
68947c6ae99SBarry Smith 
69047c6ae99SBarry Smith   /* If not Z periodic */
6911321219cSEthan Coon   if (bz != DMDA_BOUNDARY_PERIODIC) {
69247c6ae99SBarry Smith     if (zs==0)   {n0  = n1  = n2  = n3  = n4  = n5  = n6  = n7  = n8  = -2;}
69347c6ae99SBarry Smith     if (ze==P)   {n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;}
69447c6ae99SBarry Smith   }
69547c6ae99SBarry Smith 
69647c6ae99SBarry Smith   ierr = PetscMalloc(27*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr);
69747c6ae99SBarry Smith   dd->neighbors[0] = n0;
69847c6ae99SBarry Smith   dd->neighbors[1] = n1;
69947c6ae99SBarry Smith   dd->neighbors[2] = n2;
70047c6ae99SBarry Smith   dd->neighbors[3] = n3;
70147c6ae99SBarry Smith   dd->neighbors[4] = n4;
70247c6ae99SBarry Smith   dd->neighbors[5] = n5;
70347c6ae99SBarry Smith   dd->neighbors[6] = n6;
70447c6ae99SBarry Smith   dd->neighbors[7] = n7;
70547c6ae99SBarry Smith   dd->neighbors[8] = n8;
70647c6ae99SBarry Smith   dd->neighbors[9] = n9;
70747c6ae99SBarry Smith   dd->neighbors[10] = n10;
70847c6ae99SBarry Smith   dd->neighbors[11] = n11;
70947c6ae99SBarry Smith   dd->neighbors[12] = n12;
71047c6ae99SBarry Smith   dd->neighbors[13] = rank;
71147c6ae99SBarry Smith   dd->neighbors[14] = n14;
71247c6ae99SBarry Smith   dd->neighbors[15] = n15;
71347c6ae99SBarry Smith   dd->neighbors[16] = n16;
71447c6ae99SBarry Smith   dd->neighbors[17] = n17;
71547c6ae99SBarry Smith   dd->neighbors[18] = n18;
71647c6ae99SBarry Smith   dd->neighbors[19] = n19;
71747c6ae99SBarry Smith   dd->neighbors[20] = n20;
71847c6ae99SBarry Smith   dd->neighbors[21] = n21;
71947c6ae99SBarry Smith   dd->neighbors[22] = n22;
72047c6ae99SBarry Smith   dd->neighbors[23] = n23;
72147c6ae99SBarry Smith   dd->neighbors[24] = n24;
72247c6ae99SBarry Smith   dd->neighbors[25] = n25;
72347c6ae99SBarry Smith   dd->neighbors[26] = n26;
72447c6ae99SBarry Smith 
72547c6ae99SBarry Smith   /* If star stencil then delete the corner neighbors */
726aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_STAR) {
72747c6ae99SBarry Smith      /* save information about corner neighbors */
72847c6ae99SBarry Smith      sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7;
72947c6ae99SBarry Smith      sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18;
73047c6ae99SBarry Smith      sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25;
73147c6ae99SBarry Smith      sn26 = n26;
73247c6ae99SBarry Smith      n0  = n1  = n2  = n3  = n5  = n6  = n7  = n8  = n9  = n11 =
73347c6ae99SBarry Smith      n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1;
73447c6ae99SBarry Smith   }
73547c6ae99SBarry Smith 
73647c6ae99SBarry Smith 
73747c6ae99SBarry Smith   ierr = PetscMalloc((Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
73847c6ae99SBarry Smith   ierr = PetscLogObjectMemory(da,(Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(PetscInt));CHKERRQ(ierr);
73947c6ae99SBarry Smith 
74047c6ae99SBarry Smith   nn = 0;
74147c6ae99SBarry Smith   /* Bottom Level */
74247c6ae99SBarry Smith   for (k=0; k<s_z; k++) {
74347c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
74447c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
745ce00eea3SSatish Balay         x_t = lx[n0 % m];
74647c6ae99SBarry Smith         y_t = ly[(n0 % (m*n))/m];
74747c6ae99SBarry Smith         z_t = lz[n0 / (m*n)];
74847c6ae99SBarry 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;
749fdc81ce8SEthan Coon         if (s_t < 0) {s_t = bases[n0] + x_t*y_t*z_t - (s_y-i)*x_t - s_x;} /* 2D case */
75047c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
75147c6ae99SBarry Smith       }
75247c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
75347c6ae99SBarry Smith         x_t = x;
75447c6ae99SBarry Smith         y_t = ly[(n1 % (m*n))/m];
75547c6ae99SBarry Smith         z_t = lz[n1 / (m*n)];
75647c6ae99SBarry 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;
757fdc81ce8SEthan Coon         if (s_t < 0) {s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t;} /* 2D case */
75847c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
75947c6ae99SBarry Smith       }
76047c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
761ce00eea3SSatish Balay         x_t = lx[n2 % m];
76247c6ae99SBarry Smith         y_t = ly[(n2 % (m*n))/m];
76347c6ae99SBarry Smith         z_t = lz[n2 / (m*n)];
76447c6ae99SBarry 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;
765fdc81ce8SEthan Coon         if (s_t < 0) {s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t;} /* 2D case */
76647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
76747c6ae99SBarry Smith       }
76847c6ae99SBarry Smith     }
76947c6ae99SBarry Smith 
77047c6ae99SBarry Smith     for (i=0; i<y; i++) {
77147c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
772ce00eea3SSatish Balay         x_t = lx[n3 % m];
77347c6ae99SBarry Smith         y_t = y;
77447c6ae99SBarry Smith         z_t = lz[n3 / (m*n)];
77547c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
776fdc81ce8SEthan Coon         if (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 */
77747c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
77847c6ae99SBarry Smith       }
77947c6ae99SBarry Smith 
78047c6ae99SBarry Smith       if (n4 >= 0) { /* middle */
78147c6ae99SBarry Smith         x_t = x;
78247c6ae99SBarry Smith         y_t = y;
78347c6ae99SBarry Smith         z_t = lz[n4 / (m*n)];
78447c6ae99SBarry Smith         s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
785fdc81ce8SEthan Coon         if (s_t < 0) {s_t = bases[n4] + i*x_t + x_t*y_t*z_t - x_t*y_t;} /* 2D case */
78647c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
78747c6ae99SBarry Smith       }
78847c6ae99SBarry Smith 
78947c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
790ce00eea3SSatish Balay         x_t = lx[n5 % m];
79147c6ae99SBarry Smith         y_t = y;
79247c6ae99SBarry Smith         z_t = lz[n5 / (m*n)];
79347c6ae99SBarry Smith         s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
794fdc81ce8SEthan Coon         if (s_t < 0) {s_t = bases[n5] + i*x_t + x_t*y_t*z_t - x_t*y_t;} /* 2D case */
79547c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
79647c6ae99SBarry Smith       }
79747c6ae99SBarry Smith     }
79847c6ae99SBarry Smith 
79947c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
80047c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
801ce00eea3SSatish Balay         x_t = lx[n6 % m];
80247c6ae99SBarry Smith         y_t = ly[(n6 % (m*n))/m];
80347c6ae99SBarry Smith         z_t = lz[n6 / (m*n)];
80447c6ae99SBarry Smith         s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
805fdc81ce8SEthan Coon         if (s_t < 0) {s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - x_t*y_t;} /* 2D case */
80647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
80747c6ae99SBarry Smith       }
80847c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
80947c6ae99SBarry Smith         x_t = x;
81047c6ae99SBarry Smith         y_t = ly[(n7 % (m*n))/m];
81147c6ae99SBarry Smith         z_t = lz[n7 / (m*n)];
81247c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
813fdc81ce8SEthan Coon         if (s_t < 0) {s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - x_t*y_t;} /* 2D case */
81447c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
81547c6ae99SBarry Smith       }
81647c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
817ce00eea3SSatish Balay         x_t = lx[n8 % m];
81847c6ae99SBarry Smith         y_t = ly[(n8 % (m*n))/m];
81947c6ae99SBarry Smith         z_t = lz[n8 / (m*n)];
82047c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
821fdc81ce8SEthan Coon         if (s_t < 0) {s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - x_t*y_t;} /* 2D case */
82247c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
82347c6ae99SBarry Smith       }
82447c6ae99SBarry Smith     }
82547c6ae99SBarry Smith   }
82647c6ae99SBarry Smith 
82747c6ae99SBarry Smith   /* Middle Level */
82847c6ae99SBarry Smith   for (k=0; k<z; k++) {
82947c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
83047c6ae99SBarry Smith       if (n9 >= 0) { /* left below */
831ce00eea3SSatish Balay         x_t = lx[n9 % m];
83247c6ae99SBarry Smith         y_t = ly[(n9 % (m*n))/m];
83347c6ae99SBarry Smith         /* z_t = z; */
83447c6ae99SBarry Smith         s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
83547c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
83647c6ae99SBarry Smith       }
83747c6ae99SBarry Smith       if (n10 >= 0) { /* directly below */
83847c6ae99SBarry Smith         x_t = x;
83947c6ae99SBarry Smith         y_t = ly[(n10 % (m*n))/m];
84047c6ae99SBarry Smith         /* z_t = z; */
84147c6ae99SBarry Smith         s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
84247c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
84347c6ae99SBarry Smith       }
84447c6ae99SBarry Smith       if (n11 >= 0) { /* right below */
845ce00eea3SSatish Balay         x_t = lx[n11 % m];
84647c6ae99SBarry Smith         y_t = ly[(n11 % (m*n))/m];
84747c6ae99SBarry Smith         /* z_t = z; */
84847c6ae99SBarry Smith         s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
84947c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
85047c6ae99SBarry Smith       }
85147c6ae99SBarry Smith     }
85247c6ae99SBarry Smith 
85347c6ae99SBarry Smith     for (i=0; i<y; i++) {
85447c6ae99SBarry Smith       if (n12 >= 0) { /* directly left */
855ce00eea3SSatish Balay         x_t = lx[n12 % m];
85647c6ae99SBarry Smith         y_t = y;
85747c6ae99SBarry Smith         /* z_t = z; */
85847c6ae99SBarry Smith         s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
85947c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
86047c6ae99SBarry Smith       }
86147c6ae99SBarry Smith 
86247c6ae99SBarry Smith       /* Interior */
86347c6ae99SBarry Smith       s_t = bases[rank] + i*x + k*x*y;
86447c6ae99SBarry Smith       for (j=0; j<x; j++) { idx[nn++] = s_t++;}
86547c6ae99SBarry Smith 
86647c6ae99SBarry Smith       if (n14 >= 0) { /* directly right */
867ce00eea3SSatish Balay         x_t = lx[n14 % m];
86847c6ae99SBarry Smith         y_t = y;
86947c6ae99SBarry Smith         /* z_t = z; */
87047c6ae99SBarry Smith         s_t = bases[n14] + i*x_t + k*x_t*y_t;
87147c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
87247c6ae99SBarry Smith       }
87347c6ae99SBarry Smith     }
87447c6ae99SBarry Smith 
87547c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
87647c6ae99SBarry Smith       if (n15 >= 0) { /* left above */
877ce00eea3SSatish Balay         x_t = lx[n15 % m];
87847c6ae99SBarry Smith         y_t = ly[(n15 % (m*n))/m];
87947c6ae99SBarry Smith         /* z_t = z; */
88047c6ae99SBarry Smith         s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
88147c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
88247c6ae99SBarry Smith       }
88347c6ae99SBarry Smith       if (n16 >= 0) { /* directly above */
88447c6ae99SBarry Smith         x_t = x;
88547c6ae99SBarry Smith         y_t = ly[(n16 % (m*n))/m];
88647c6ae99SBarry Smith         /* z_t = z; */
88747c6ae99SBarry Smith         s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
88847c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
88947c6ae99SBarry Smith       }
89047c6ae99SBarry Smith       if (n17 >= 0) { /* right above */
891ce00eea3SSatish Balay         x_t = lx[n17 % m];
89247c6ae99SBarry Smith         y_t = ly[(n17 % (m*n))/m];
89347c6ae99SBarry Smith         /* z_t = z; */
89447c6ae99SBarry Smith         s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
89547c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
89647c6ae99SBarry Smith       }
89747c6ae99SBarry Smith     }
89847c6ae99SBarry Smith   }
89947c6ae99SBarry Smith 
90047c6ae99SBarry Smith   /* Upper Level */
90147c6ae99SBarry Smith   for (k=0; k<s_z; k++) {
90247c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
90347c6ae99SBarry Smith       if (n18 >= 0) { /* left below */
904ce00eea3SSatish Balay         x_t = lx[n18 % m];
90547c6ae99SBarry Smith         y_t = ly[(n18 % (m*n))/m];
90647c6ae99SBarry Smith         /* z_t = lz[n18 / (m*n)]; */
90747c6ae99SBarry Smith         s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
908fdc81ce8SEthan Coon         if (s_t >= x*y*z) {s_t = bases[n18] - (s_y-i)*x_t -s_x + x_t*y_t;} /* 2d case */
90947c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
91047c6ae99SBarry Smith       }
91147c6ae99SBarry Smith       if (n19 >= 0) { /* directly below */
91247c6ae99SBarry Smith         x_t = x;
91347c6ae99SBarry Smith         y_t = ly[(n19 % (m*n))/m];
91447c6ae99SBarry Smith         /* z_t = lz[n19 / (m*n)]; */
91547c6ae99SBarry Smith         s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
916fdc81ce8SEthan Coon         if (s_t >= x*y*z) {s_t = bases[n19] - (s_y+1-i)*x_t + x_t*y_t;} /* 2d case */
91747c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
91847c6ae99SBarry Smith       }
91947c6ae99SBarry Smith       if (n20 >= 0) { /* right below */
920ce00eea3SSatish Balay         x_t = lx[n20 % m];
92147c6ae99SBarry Smith         y_t = ly[(n20 % (m*n))/m];
92247c6ae99SBarry Smith         /* z_t = lz[n20 / (m*n)]; */
92347c6ae99SBarry Smith         s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
924fdc81ce8SEthan Coon         if (s_t >= x*y*z) {s_t = bases[n20] - (s_y+1-i)*x_t + x_t*y_t;} /* 2d case */
92547c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
92647c6ae99SBarry Smith       }
92747c6ae99SBarry Smith     }
92847c6ae99SBarry Smith 
92947c6ae99SBarry Smith     for (i=0; i<y; i++) {
93047c6ae99SBarry Smith       if (n21 >= 0) { /* directly left */
931ce00eea3SSatish Balay         x_t = lx[n21 % m];
93247c6ae99SBarry Smith         y_t = y;
93347c6ae99SBarry Smith         /* z_t = lz[n21 / (m*n)]; */
93447c6ae99SBarry Smith         s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
935fdc81ce8SEthan Coon         if (s_t >= x*y*z) {s_t = bases[n21] + (i+1)*x_t - s_x;}  /* 2d case */
93647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
93747c6ae99SBarry Smith       }
93847c6ae99SBarry Smith 
93947c6ae99SBarry Smith       if (n22 >= 0) { /* middle */
94047c6ae99SBarry Smith         x_t = x;
94147c6ae99SBarry Smith         y_t = y;
94247c6ae99SBarry Smith         /* z_t = lz[n22 / (m*n)]; */
94347c6ae99SBarry Smith         s_t = bases[n22] + i*x_t + k*x_t*y_t;
944fdc81ce8SEthan Coon         if (s_t >= x*y*z) {s_t = bases[n22] + i*x_t;} /* 2d case */
94547c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
94647c6ae99SBarry Smith       }
94747c6ae99SBarry Smith 
94847c6ae99SBarry Smith       if (n23 >= 0) { /* directly right */
949ce00eea3SSatish Balay         x_t = lx[n23 % m];
95047c6ae99SBarry Smith         y_t = y;
95147c6ae99SBarry Smith         /* z_t = lz[n23 / (m*n)]; */
95247c6ae99SBarry Smith         s_t = bases[n23] + i*x_t + k*x_t*y_t;
953fdc81ce8SEthan Coon         if (s_t >= x*y*z) {s_t = bases[n23] + i*x_t;} /* 2d case */
95447c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
95547c6ae99SBarry Smith       }
95647c6ae99SBarry Smith     }
95747c6ae99SBarry Smith 
95847c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
95947c6ae99SBarry Smith       if (n24 >= 0) { /* left above */
960ce00eea3SSatish Balay         x_t = lx[n24 % m];
96147c6ae99SBarry Smith         y_t = ly[(n24 % (m*n))/m];
96247c6ae99SBarry Smith         /* z_t = lz[n24 / (m*n)]; */
96347c6ae99SBarry Smith         s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
964fdc81ce8SEthan Coon         if (s_t >= x*y*z) {s_t = bases[n24] + i*x_t - s_x;} /* 2d case */
96547c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
96647c6ae99SBarry Smith       }
96747c6ae99SBarry Smith       if (n25 >= 0) { /* directly above */
96847c6ae99SBarry Smith         x_t = x;
96947c6ae99SBarry Smith         y_t = ly[(n25 % (m*n))/m];
97047c6ae99SBarry Smith         /* z_t = lz[n25 / (m*n)]; */
97147c6ae99SBarry Smith         s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
972fdc81ce8SEthan Coon         if (s_t >= x*y*z) {s_t = bases[n25] + (i-1)*x_t;} /* 2d case */
97347c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
97447c6ae99SBarry Smith       }
97547c6ae99SBarry Smith       if (n26 >= 0) { /* right above */
976ce00eea3SSatish Balay         x_t = lx[n26 % m];
97747c6ae99SBarry Smith         y_t = ly[(n26 % (m*n))/m];
97847c6ae99SBarry Smith         /* z_t = lz[n26 / (m*n)]; */
97947c6ae99SBarry Smith         s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
980fdc81ce8SEthan Coon         if (s_t >= x*y*z) {s_t = bases[n26] + (i-1)*x_t;} /* 2d case */
98147c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
98247c6ae99SBarry Smith       }
98347c6ae99SBarry Smith     }
98447c6ae99SBarry Smith   }
985ce00eea3SSatish Balay 
986ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
98747c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
98847c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr);
989fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
990fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
99147c6ae99SBarry Smith 
992aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_STAR) {
99347c6ae99SBarry Smith     n0  = sn0;  n1  = sn1;  n2  = sn2;  n3  = sn3;  n5  = sn5;  n6  = sn6; n7 = sn7;
99447c6ae99SBarry Smith     n8  = sn8;  n9  = sn9;  n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18;
99547c6ae99SBarry Smith     n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25;
99647c6ae99SBarry Smith     n26 = sn26;
997ce00eea3SSatish Balay   }
99847c6ae99SBarry Smith 
999ce00eea3SSatish Balay   if ((stencil_type == DMDA_STENCIL_STAR) ||
10001321219cSEthan Coon       (bx != DMDA_BOUNDARY_PERIODIC && bx) ||
10011321219cSEthan Coon       (by != DMDA_BOUNDARY_PERIODIC && by) ||
10021321219cSEthan Coon       (bz != DMDA_BOUNDARY_PERIODIC && bz)) {
1003ce00eea3SSatish Balay     /*
1004ce00eea3SSatish Balay         Recompute the local to global mappings, this time keeping the
1005ce00eea3SSatish Balay       information about the cross corner processor numbers.
1006ce00eea3SSatish Balay     */
100747c6ae99SBarry Smith     nn = 0;
100847c6ae99SBarry Smith     /* Bottom Level */
100947c6ae99SBarry Smith     for (k=0; k<s_z; k++) {
101047c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
101147c6ae99SBarry Smith         if (n0 >= 0) { /* left below */
1012ce00eea3SSatish Balay           x_t = lx[n0 % m];
101347c6ae99SBarry Smith           y_t = ly[(n0 % (m*n))/m];
101447c6ae99SBarry Smith           z_t = lz[n0 / (m*n)];
101547c6ae99SBarry 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;
101647c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1017ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Ys-ys < 0 && Zs-zs < 0) {
1018ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
101947c6ae99SBarry Smith         }
102047c6ae99SBarry Smith         if (n1 >= 0) { /* directly below */
102147c6ae99SBarry Smith           x_t = x;
102247c6ae99SBarry Smith           y_t = ly[(n1 % (m*n))/m];
102347c6ae99SBarry Smith           z_t = lz[n1 / (m*n)];
102447c6ae99SBarry 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;
102547c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1026ce00eea3SSatish Balay         } else if (Ys-ys < 0 && Zs-zs < 0) {
1027ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
102847c6ae99SBarry Smith         }
102947c6ae99SBarry Smith         if (n2 >= 0) { /* right below */
1030ce00eea3SSatish Balay           x_t = lx[n2 % m];
103147c6ae99SBarry Smith           y_t = ly[(n2 % (m*n))/m];
103247c6ae99SBarry Smith           z_t = lz[n2 / (m*n)];
103347c6ae99SBarry 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;
103447c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1035ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Ys-ys < 0 && Zs-zs < 0) {
1036ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
103747c6ae99SBarry Smith         }
103847c6ae99SBarry Smith       }
103947c6ae99SBarry Smith 
104047c6ae99SBarry Smith       for (i=0; i<y; i++) {
104147c6ae99SBarry Smith         if (n3 >= 0) { /* directly left */
1042ce00eea3SSatish Balay           x_t = lx[n3 % m];
104347c6ae99SBarry Smith           y_t = y;
104447c6ae99SBarry Smith           z_t = lz[n3 / (m*n)];
104547c6ae99SBarry Smith           s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
104647c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1047ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Zs-zs < 0) {
1048ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
104947c6ae99SBarry Smith         }
105047c6ae99SBarry Smith 
105147c6ae99SBarry Smith         if (n4 >= 0) { /* middle */
105247c6ae99SBarry Smith           x_t = x;
105347c6ae99SBarry Smith           y_t = y;
105447c6ae99SBarry Smith           z_t = lz[n4 / (m*n)];
105547c6ae99SBarry Smith           s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
105647c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1057ce00eea3SSatish Balay         } else if (Zs-zs < 0) {
1058ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
105947c6ae99SBarry Smith         }
106047c6ae99SBarry Smith 
106147c6ae99SBarry Smith         if (n5 >= 0) { /* directly right */
1062ce00eea3SSatish Balay           x_t = lx[n5 % m];
106347c6ae99SBarry Smith           y_t = y;
106447c6ae99SBarry Smith           z_t = lz[n5 / (m*n)];
106547c6ae99SBarry Smith           s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
106647c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1067ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Zs-zs < 0) {
1068ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
106947c6ae99SBarry Smith         }
107047c6ae99SBarry Smith       }
107147c6ae99SBarry Smith 
107247c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
107347c6ae99SBarry Smith         if (n6 >= 0) { /* left above */
1074ce00eea3SSatish Balay           x_t = lx[n6 % m];
107547c6ae99SBarry Smith           y_t = ly[(n6 % (m*n))/m];
107647c6ae99SBarry Smith           z_t = lz[n6 / (m*n)];
107747c6ae99SBarry Smith           s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
107847c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1079ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ye-Ye < 0 && Zs-zs < 0) {
1080ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
108147c6ae99SBarry Smith         }
108247c6ae99SBarry Smith         if (n7 >= 0) { /* directly above */
108347c6ae99SBarry Smith           x_t = x;
108447c6ae99SBarry Smith           y_t = ly[(n7 % (m*n))/m];
108547c6ae99SBarry Smith           z_t = lz[n7 / (m*n)];
108647c6ae99SBarry Smith           s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
108747c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1088ce00eea3SSatish Balay         } else if (ye-Ye < 0 && Zs-zs < 0) {
1089ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
109047c6ae99SBarry Smith         }
109147c6ae99SBarry Smith         if (n8 >= 0) { /* right above */
1092ce00eea3SSatish Balay           x_t = lx[n8 % m];
109347c6ae99SBarry Smith           y_t = ly[(n8 % (m*n))/m];
109447c6ae99SBarry Smith           z_t = lz[n8 / (m*n)];
109547c6ae99SBarry Smith           s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
109647c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1097ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ye-Ye < 0 && Zs-zs < 0) {
1098ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
109947c6ae99SBarry Smith         }
110047c6ae99SBarry Smith       }
110147c6ae99SBarry Smith     }
110247c6ae99SBarry Smith 
110347c6ae99SBarry Smith     /* Middle Level */
110447c6ae99SBarry Smith     for (k=0; k<z; k++) {
110547c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
110647c6ae99SBarry Smith         if (n9 >= 0) { /* left below */
1107ce00eea3SSatish Balay           x_t = lx[n9 % m];
110847c6ae99SBarry Smith           y_t = ly[(n9 % (m*n))/m];
110947c6ae99SBarry Smith           /* z_t = z; */
111047c6ae99SBarry Smith           s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
111147c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1112ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Ys-ys < 0) {
1113ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
111447c6ae99SBarry Smith         }
111547c6ae99SBarry Smith         if (n10 >= 0) { /* directly below */
111647c6ae99SBarry Smith           x_t = x;
111747c6ae99SBarry Smith           y_t = ly[(n10 % (m*n))/m];
111847c6ae99SBarry Smith           /* z_t = z; */
111947c6ae99SBarry Smith           s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
112047c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1121ce00eea3SSatish Balay         } else if (Ys-ys < 0) {
1122ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
112347c6ae99SBarry Smith         }
112447c6ae99SBarry Smith         if (n11 >= 0) { /* right below */
1125ce00eea3SSatish Balay           x_t = lx[n11 % m];
112647c6ae99SBarry Smith           y_t = ly[(n11 % (m*n))/m];
112747c6ae99SBarry Smith           /* z_t = z; */
112847c6ae99SBarry Smith           s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
112947c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1130ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Ys-ys < 0) {
1131ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
113247c6ae99SBarry Smith         }
113347c6ae99SBarry Smith       }
113447c6ae99SBarry Smith 
113547c6ae99SBarry Smith       for (i=0; i<y; i++) {
113647c6ae99SBarry Smith         if (n12 >= 0) { /* directly left */
1137ce00eea3SSatish Balay           x_t = lx[n12 % m];
113847c6ae99SBarry Smith           y_t = y;
113947c6ae99SBarry Smith           /* z_t = z; */
114047c6ae99SBarry Smith           s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
114147c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1142ce00eea3SSatish Balay         } else if (Xs-xs < 0) {
1143ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
114447c6ae99SBarry Smith         }
114547c6ae99SBarry Smith 
114647c6ae99SBarry Smith         /* Interior */
114747c6ae99SBarry Smith         s_t = bases[rank] + i*x + k*x*y;
114847c6ae99SBarry Smith         for (j=0; j<x; j++) { idx[nn++] = s_t++;}
114947c6ae99SBarry Smith 
115047c6ae99SBarry Smith         if (n14 >= 0) { /* directly right */
1151ce00eea3SSatish Balay           x_t = lx[n14 % m];
115247c6ae99SBarry Smith           y_t = y;
115347c6ae99SBarry Smith           /* z_t = z; */
115447c6ae99SBarry Smith           s_t = bases[n14] + i*x_t + k*x_t*y_t;
115547c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1156ce00eea3SSatish Balay         } else if (xe-Xe < 0) {
1157ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
115847c6ae99SBarry Smith         }
115947c6ae99SBarry Smith       }
116047c6ae99SBarry Smith 
116147c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
116247c6ae99SBarry Smith         if (n15 >= 0) { /* left above */
1163ce00eea3SSatish Balay           x_t = lx[n15 % m];
116447c6ae99SBarry Smith           y_t = ly[(n15 % (m*n))/m];
116547c6ae99SBarry Smith           /* z_t = z; */
116647c6ae99SBarry Smith           s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
116747c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1168ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ye-Ye < 0) {
1169ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
117047c6ae99SBarry Smith         }
117147c6ae99SBarry Smith         if (n16 >= 0) { /* directly above */
117247c6ae99SBarry Smith           x_t = x;
117347c6ae99SBarry Smith           y_t = ly[(n16 % (m*n))/m];
117447c6ae99SBarry Smith           /* z_t = z; */
117547c6ae99SBarry Smith           s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
117647c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1177ce00eea3SSatish Balay         } else if (ye-Ye < 0) {
1178ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
117947c6ae99SBarry Smith         }
118047c6ae99SBarry Smith         if (n17 >= 0) { /* right above */
1181ce00eea3SSatish Balay           x_t = lx[n17 % m];
118247c6ae99SBarry Smith           y_t = ly[(n17 % (m*n))/m];
118347c6ae99SBarry Smith           /* z_t = z; */
118447c6ae99SBarry Smith           s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
118547c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1186ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ye-Ye < 0) {
1187ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
118847c6ae99SBarry Smith         }
118947c6ae99SBarry Smith       }
119047c6ae99SBarry Smith     }
119147c6ae99SBarry Smith 
119247c6ae99SBarry Smith     /* Upper Level */
119347c6ae99SBarry Smith     for (k=0; k<s_z; k++) {
119447c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
119547c6ae99SBarry Smith         if (n18 >= 0) { /* left below */
1196ce00eea3SSatish Balay           x_t = lx[n18 % m];
119747c6ae99SBarry Smith           y_t = ly[(n18 % (m*n))/m];
119847c6ae99SBarry Smith           /* z_t = lz[n18 / (m*n)]; */
119947c6ae99SBarry Smith           s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
120047c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1201ce00eea3SSatish Balay         } else if (Xs-xs < 0 && Ys-ys < 0 && ze-Ze < 0) {
1202ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
120347c6ae99SBarry Smith         }
120447c6ae99SBarry Smith         if (n19 >= 0) { /* directly below */
120547c6ae99SBarry Smith           x_t = x;
120647c6ae99SBarry Smith           y_t = ly[(n19 % (m*n))/m];
120747c6ae99SBarry Smith           /* z_t = lz[n19 / (m*n)]; */
120847c6ae99SBarry Smith           s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
120947c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1210ce00eea3SSatish Balay         } else if (Ys-ys < 0 && ze-Ze < 0) {
1211ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
121247c6ae99SBarry Smith         }
121347c6ae99SBarry Smith         if (n20 >= 0) { /* right below */
1214ce00eea3SSatish Balay           x_t = lx[n20 % m];
121547c6ae99SBarry Smith           y_t = ly[(n20 % (m*n))/m];
121647c6ae99SBarry Smith           /* z_t = lz[n20 / (m*n)]; */
121747c6ae99SBarry Smith           s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
121847c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1219ce00eea3SSatish Balay         } else if (xe-Xe < 0 && Ys-ys < 0 && ze-Ze < 0) {
1220ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
122147c6ae99SBarry Smith         }
122247c6ae99SBarry Smith       }
122347c6ae99SBarry Smith 
122447c6ae99SBarry Smith       for (i=0; i<y; i++) {
122547c6ae99SBarry Smith         if (n21 >= 0) { /* directly left */
1226ce00eea3SSatish Balay           x_t = lx[n21 % m];
122747c6ae99SBarry Smith           y_t = y;
122847c6ae99SBarry Smith           /* z_t = lz[n21 / (m*n)]; */
122947c6ae99SBarry Smith           s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
123047c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1231ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ze-Ze < 0) {
1232ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
123347c6ae99SBarry Smith         }
123447c6ae99SBarry Smith 
123547c6ae99SBarry Smith         if (n22 >= 0) { /* middle */
123647c6ae99SBarry Smith           x_t = x;
123747c6ae99SBarry Smith           y_t = y;
123847c6ae99SBarry Smith           /* z_t = lz[n22 / (m*n)]; */
123947c6ae99SBarry Smith           s_t = bases[n22] + i*x_t + k*x_t*y_t;
124047c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1241ce00eea3SSatish Balay         } else if (ze-Ze < 0) {
1242ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
124347c6ae99SBarry Smith         }
124447c6ae99SBarry Smith 
124547c6ae99SBarry Smith         if (n23 >= 0) { /* directly right */
1246ce00eea3SSatish Balay           x_t = lx[n23 % m];
124747c6ae99SBarry Smith           y_t = y;
124847c6ae99SBarry Smith           /* z_t = lz[n23 / (m*n)]; */
124947c6ae99SBarry Smith           s_t = bases[n23] + i*x_t + k*x_t*y_t;
125047c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1251ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ze-Ze < 0) {
1252ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
125347c6ae99SBarry Smith         }
125447c6ae99SBarry Smith       }
125547c6ae99SBarry Smith 
125647c6ae99SBarry Smith       for (i=1; i<=s_y; i++) {
125747c6ae99SBarry Smith         if (n24 >= 0) { /* left above */
1258ce00eea3SSatish Balay           x_t = lx[n24 % m];
125947c6ae99SBarry Smith           y_t = ly[(n24 % (m*n))/m];
126047c6ae99SBarry Smith           /* z_t = lz[n24 / (m*n)]; */
126147c6ae99SBarry Smith           s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
126247c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1263ce00eea3SSatish Balay         } else if (Xs-xs < 0 && ye-Ye < 0 && ze-Ze < 0) {
1264ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
126547c6ae99SBarry Smith         }
126647c6ae99SBarry Smith         if (n25 >= 0) { /* directly above */
126747c6ae99SBarry Smith           x_t = x;
126847c6ae99SBarry Smith           y_t = ly[(n25 % (m*n))/m];
126947c6ae99SBarry Smith           /* z_t = lz[n25 / (m*n)]; */
127047c6ae99SBarry Smith           s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
127147c6ae99SBarry Smith           for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1272ce00eea3SSatish Balay         } else if (ye-Ye < 0 && ze-Ze < 0) {
1273ce00eea3SSatish Balay           for (j=0; j<x; j++) { idx[nn++] = -1;}
127447c6ae99SBarry Smith         }
127547c6ae99SBarry Smith         if (n26 >= 0) { /* right above */
1276ce00eea3SSatish Balay           x_t = lx[n26 % m];
127747c6ae99SBarry Smith           y_t = ly[(n26 % (m*n))/m];
127847c6ae99SBarry Smith           /* z_t = lz[n26 / (m*n)]; */
127947c6ae99SBarry Smith           s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
128047c6ae99SBarry Smith           for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1281ce00eea3SSatish Balay         } else if (xe-Xe < 0 && ye-Ye < 0 && ze-Ze < 0) {
1282ce00eea3SSatish Balay           for (j=0; j<s_x; j++) { idx[nn++] = -1;}
128347c6ae99SBarry Smith         }
128447c6ae99SBarry Smith       }
128547c6ae99SBarry Smith     }
128647c6ae99SBarry Smith   }
128747c6ae99SBarry Smith   /*
128847c6ae99SBarry Smith      Set the local to global ordering in the global vector, this allows use
128947c6ae99SBarry Smith      of VecSetValuesLocal().
129047c6ae99SBarry Smith   */
1291ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_OWN_POINTER,&ltogis);CHKERRQ(ierr);
1292ce00eea3SSatish Balay   ierr = PetscMalloc(nn*dof*sizeof(PetscInt),&idx_cpy);CHKERRQ(ierr);
1293db87c5ecSEthan Coon   ierr = PetscLogObjectMemory(da,nn*dof*sizeof(PetscInt));CHKERRQ(ierr);
1294ce00eea3SSatish Balay   ierr = ISGetIndices(ltogis, &idx_full);
1295ce00eea3SSatish Balay   ierr = PetscMemcpy(idx_cpy,idx_full,nn*dof*sizeof(PetscInt));CHKERRQ(ierr);
1296ce00eea3SSatish Balay   ierr = ISRestoreIndices(ltogis, &idx_full);
1297ce00eea3SSatish Balay   ierr = ISLocalToGlobalMappingCreateIS(ltogis,&da->ltogmap);CHKERRQ(ierr);
1298ce00eea3SSatish Balay   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
1299fcfd50ebSBarry Smith   ierr = ISDestroy(&ltogis);CHKERRQ(ierr);
13001411c6eeSJed Brown   ierr = ISLocalToGlobalMappingBlock(da->ltogmap,dd->w,&da->ltogmapb);CHKERRQ(ierr);
13011411c6eeSJed Brown   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
130247c6ae99SBarry Smith 
1303ce00eea3SSatish Balay   ierr = PetscFree2(bases,ldims);CHKERRQ(ierr);
1304ce00eea3SSatish Balay   dd->m  = m;  dd->n  = n;  dd->p  = p;
1305ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
1306ce00eea3SSatish Balay   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = zs; dd->ze = ze;
1307ce00eea3SSatish Balay   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = Zs; dd->Ze = Ze;
1308ce00eea3SSatish Balay 
1309fcfd50ebSBarry Smith   ierr = VecDestroy(&local);CHKERRQ(ierr);
1310fcfd50ebSBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
1311ce00eea3SSatish Balay 
1312ce00eea3SSatish Balay   dd->gtol      = gtol;
1313ce00eea3SSatish Balay   dd->ltog      = ltog;
1314ce00eea3SSatish Balay   dd->idx       = idx_cpy;
1315ce00eea3SSatish Balay   dd->Nl        = nn*dof;
1316ce00eea3SSatish Balay   dd->base      = base;
1317ce00eea3SSatish Balay   da->ops->view = DMView_DA_3d;
131847c6ae99SBarry Smith   dd->ltol = PETSC_NULL;
131947c6ae99SBarry Smith   dd->ao   = PETSC_NULL;
1320ce00eea3SSatish Balay 
132147c6ae99SBarry Smith   PetscFunctionReturn(0);
132247c6ae99SBarry Smith }
1323564755cdSBarry Smith 
132447c6ae99SBarry Smith 
132547c6ae99SBarry Smith #undef __FUNCT__
1326aa219208SBarry Smith #define __FUNCT__ "DMDACreate3d"
132747c6ae99SBarry Smith /*@C
1328aa219208SBarry Smith    DMDACreate3d - Creates an object that will manage the communication of three-dimensional
132947c6ae99SBarry Smith    regular array data that is distributed across some processors.
133047c6ae99SBarry Smith 
133147c6ae99SBarry Smith    Collective on MPI_Comm
133247c6ae99SBarry Smith 
133347c6ae99SBarry Smith    Input Parameters:
133447c6ae99SBarry Smith +  comm - MPI communicator
13351321219cSEthan Coon .  bx,by,bz - type of ghost nodes the array have.
13361321219cSEthan Coon          Use one of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC.
1337aa219208SBarry Smith .  stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX)
133847c6ae99SBarry 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
133947c6ae99SBarry Smith             from the command line with -da_grid_x <M> -da_grid_y <N> -da_grid_z <P>)
134047c6ae99SBarry Smith .  m,n,p - corresponding number of processors in each dimension
134147c6ae99SBarry Smith            (or PETSC_DECIDE to have calculated)
134247c6ae99SBarry Smith .  dof - number of degrees of freedom per node
134310d7c030SMatthew G Knepley .  s - stencil width
134410d7c030SMatthew G Knepley -  lx, ly, lz - arrays containing the number of nodes in each cell along
134547c6ae99SBarry Smith           the x, y, and z coordinates, or PETSC_NULL. If non-null, these
134647c6ae99SBarry Smith           must be of length as m,n,p and the corresponding
134747c6ae99SBarry Smith           m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of
134847c6ae99SBarry Smith           the ly[] must N, sum of the lz[] must be P
134947c6ae99SBarry Smith 
135047c6ae99SBarry Smith    Output Parameter:
135147c6ae99SBarry Smith .  da - the resulting distributed array object
135247c6ae99SBarry Smith 
135347c6ae99SBarry Smith    Options Database Key:
1354aa219208SBarry Smith +  -da_view - Calls DMView() at the conclusion of DMDACreate3d()
135547c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
135647c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
135747c6ae99SBarry Smith .  -da_grid_z <nz> - number of grid points in z direction, if P < 0
135847c6ae99SBarry Smith .  -da_processors_x <MX> - number of processors in x direction
135947c6ae99SBarry Smith .  -da_processors_y <MY> - number of processors in y direction
136047c6ae99SBarry Smith .  -da_processors_z <MZ> - number of processors in z direction
1361e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement ratio in x direction
1362e0f5d30fSBarry Smith .  -da_refine_y <ry> - refinement ratio in y direction
1363e0f5d30fSBarry Smith .  -da_refine_z <rz>- refinement ratio in z directio
1364e0f5d30fSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating it, , if M, N, or P < 0
136547c6ae99SBarry Smith 
136647c6ae99SBarry Smith    Level: beginner
136747c6ae99SBarry Smith 
136847c6ae99SBarry Smith    Notes:
1369aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
1370aa219208SBarry Smith    standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
137147c6ae99SBarry Smith    the standard 27-pt stencil.
137247c6ae99SBarry Smith 
1373aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
1374564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
1375564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
137647c6ae99SBarry Smith 
137747c6ae99SBarry Smith .keywords: distributed array, create, three-dimensional
137847c6ae99SBarry Smith 
1379aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate2d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
1380aa219208SBarry Smith           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDASetRefinementFactor(),
1381d461ba97SBarry Smith           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
138247c6ae99SBarry Smith 
138347c6ae99SBarry Smith @*/
13841321219cSEthan Coon PetscErrorCode  DMDACreate3d(MPI_Comm comm,DMDABoundaryType bx,DMDABoundaryType by,DMDABoundaryType bz,DMDAStencilType stencil_type,PetscInt M,
13859a42bb27SBarry 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)
138647c6ae99SBarry Smith {
138747c6ae99SBarry Smith   PetscErrorCode ierr;
138847c6ae99SBarry Smith 
138947c6ae99SBarry Smith   PetscFunctionBegin;
1390aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
1391aa219208SBarry Smith   ierr = DMDASetDim(*da, 3);CHKERRQ(ierr);
1392aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, P);CHKERRQ(ierr);
1393aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, p);CHKERRQ(ierr);
13941321219cSEthan Coon   ierr = DMDASetBoundaryType(*da, bx, by, bz);CHKERRQ(ierr);
1395aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
1396aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
1397aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
1398aa219208SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, lz);CHKERRQ(ierr);
139947c6ae99SBarry Smith   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
14009a42bb27SBarry Smith   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
14019a42bb27SBarry Smith   ierr = DMSetUp(*da);CHKERRQ(ierr);
14027242296bSJed Brown   ierr = DMView_DA_Private(*da);CHKERRQ(ierr);
140347c6ae99SBarry Smith   PetscFunctionReturn(0);
140447c6ae99SBarry Smith }
1405