xref: /petsc/src/dm/impls/da/da2.c (revision ca266f367873f97695dfce72da4e09458ac9ba8a)
19a42bb27SBarry Smith 
2b45d2f2cSJed Brown #include <petsc-private/daimpl.h>    /*I   "petscdmda.h"   I*/
347c6ae99SBarry Smith 
447c6ae99SBarry Smith #undef __FUNCT__
59a42bb27SBarry Smith #define __FUNCT__ "DMView_DA_2d"
69a42bb27SBarry Smith PetscErrorCode DMView_DA_2d(DM da,PetscViewer viewer)
747c6ae99SBarry Smith {
847c6ae99SBarry Smith   PetscErrorCode ierr;
947c6ae99SBarry Smith   PetscMPIInt    rank;
109a42bb27SBarry Smith   PetscBool      iascii,isdraw,isbinary;
1147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
129a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
139a42bb27SBarry Smith   PetscBool      ismatlab;
149a42bb27SBarry Smith #endif
1547c6ae99SBarry Smith 
1647c6ae99SBarry Smith   PetscFunctionBegin;
1747c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)da)->comm,&rank);CHKERRQ(ierr);
1847c6ae99SBarry Smith 
19251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
20251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
21251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
229a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
23251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr);
249a42bb27SBarry Smith #endif
2547c6ae99SBarry Smith   if (iascii) {
2647c6ae99SBarry Smith     PetscViewerFormat format;
2747c6ae99SBarry Smith 
2847c6ae99SBarry Smith     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
2947c6ae99SBarry Smith     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) {
30aa219208SBarry Smith       DMDALocalInfo info;
31aa219208SBarry Smith       ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
327b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
3347c6ae99SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %D N %D m %D n %D w %D s %D\n",rank,dd->M,dd->N,dd->m,dd->n,dd->w,dd->s);CHKERRQ(ierr);
3447c6ae99SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D, Y range of indices: %D %D\n",info.xs,info.xs+info.xm,info.ys,info.ys+info.ym);CHKERRQ(ierr);
3547c6ae99SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
367b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
373da9ae13SJed Brown     } else {
383da9ae13SJed Brown       ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr);
3947c6ae99SBarry Smith     }
4047c6ae99SBarry Smith   } else if (isdraw) {
4147c6ae99SBarry Smith     PetscDraw  draw;
4247c6ae99SBarry Smith     double     ymin = -1*dd->s-1,ymax = dd->N+dd->s;
4347c6ae99SBarry Smith     double     xmin = -1*dd->s-1,xmax = dd->M+dd->s;
4447c6ae99SBarry Smith     double     x,y;
4547c6ae99SBarry Smith     PetscInt   base,*idx;
4647c6ae99SBarry Smith     char       node[10];
4747c6ae99SBarry Smith     PetscBool  isnull;
4847c6ae99SBarry Smith 
4947c6ae99SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
5047c6ae99SBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
516636e97aSMatthew G Knepley     if (!da->coordinates) {
5247c6ae99SBarry Smith       ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
5347c6ae99SBarry Smith     }
5447c6ae99SBarry Smith     ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
5547c6ae99SBarry Smith 
5647c6ae99SBarry Smith     /* first processor draw all node lines */
5747c6ae99SBarry Smith     if (!rank) {
5847c6ae99SBarry Smith       ymin = 0.0; ymax = dd->N - 1;
5947c6ae99SBarry Smith       for (xmin=0; xmin<dd->M; xmin++) {
6047c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
6147c6ae99SBarry Smith       }
6247c6ae99SBarry Smith       xmin = 0.0; xmax = dd->M - 1;
6347c6ae99SBarry Smith       for (ymin=0; ymin<dd->N; ymin++) {
6447c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
6547c6ae99SBarry Smith       }
6647c6ae99SBarry Smith     }
6747c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
6847c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
6947c6ae99SBarry Smith 
7047c6ae99SBarry Smith     /* draw my box */
7147c6ae99SBarry Smith     ymin = dd->ys; ymax = dd->ye - 1; xmin = dd->xs/dd->w;
7247c6ae99SBarry Smith     xmax =(dd->xe-1)/dd->w;
7347c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
7447c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
7547c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
7647c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
7747c6ae99SBarry Smith 
7847c6ae99SBarry Smith     /* put in numbers */
7947c6ae99SBarry Smith     base = (dd->base)/dd->w;
8047c6ae99SBarry Smith     for (y=ymin; y<=ymax; y++) {
8147c6ae99SBarry Smith       for (x=xmin; x<=xmax; x++) {
8247c6ae99SBarry Smith         sprintf(node,"%d",(int)base++);
8347c6ae99SBarry Smith         ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
8447c6ae99SBarry Smith       }
8547c6ae99SBarry Smith     }
8647c6ae99SBarry Smith 
8747c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
8847c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
8947c6ae99SBarry Smith     /* overlay ghost numbers, useful for error checking */
9047c6ae99SBarry Smith     /* put in numbers */
9147c6ae99SBarry Smith 
9247c6ae99SBarry Smith     base = 0; idx = dd->idx;
9347c6ae99SBarry Smith     ymin = dd->Ys; ymax = dd->Ye; xmin = dd->Xs; xmax = dd->Xe;
9447c6ae99SBarry Smith     for (y=ymin; y<ymax; y++) {
9547c6ae99SBarry Smith       for (x=xmin; x<xmax; x++) {
9647c6ae99SBarry Smith         if ((base % dd->w) == 0) {
9747c6ae99SBarry Smith           sprintf(node,"%d",(int)(idx[base]/dd->w));
9847c6ae99SBarry Smith           ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
9947c6ae99SBarry Smith         }
10047c6ae99SBarry Smith         base++;
10147c6ae99SBarry Smith       }
10247c6ae99SBarry Smith     }
10347c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
10447c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
1059a42bb27SBarry Smith   } else if (isbinary){
1069a42bb27SBarry Smith     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
1079a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1089a42bb27SBarry Smith   } else if (ismatlab) {
1099a42bb27SBarry Smith     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
1109a42bb27SBarry Smith #endif
111aa219208SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for DMDA 1d",((PetscObject)viewer)->type_name);
11247c6ae99SBarry Smith   PetscFunctionReturn(0);
11347c6ae99SBarry Smith }
11447c6ae99SBarry Smith 
11547c6ae99SBarry Smith /*
11647c6ae99SBarry Smith       M is number of grid points
11747c6ae99SBarry Smith       m is number of processors
11847c6ae99SBarry Smith 
11947c6ae99SBarry Smith */
12047c6ae99SBarry Smith #undef __FUNCT__
121aa219208SBarry Smith #define __FUNCT__ "DMDASplitComm2d"
1227087cfbeSBarry Smith PetscErrorCode  DMDASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm)
12347c6ae99SBarry Smith {
12447c6ae99SBarry Smith   PetscErrorCode ierr;
12547c6ae99SBarry Smith   PetscInt       m,n = 0,x = 0,y = 0;
12647c6ae99SBarry Smith   PetscMPIInt    size,csize,rank;
12747c6ae99SBarry Smith 
12847c6ae99SBarry Smith   PetscFunctionBegin;
12947c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
13047c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13147c6ae99SBarry Smith 
13247c6ae99SBarry Smith   csize = 4*size;
13347c6ae99SBarry Smith   do {
13447c6ae99SBarry Smith     if (csize % 4) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot split communicator of size %d tried %d %D %D",size,csize,x,y);
13547c6ae99SBarry Smith     csize   = csize/4;
13647c6ae99SBarry Smith 
13747c6ae99SBarry Smith     m = (PetscInt)(0.5 + sqrt(((double)M)*((double)csize)/((double)N)));
13847c6ae99SBarry Smith     if (!m) m = 1;
13947c6ae99SBarry Smith     while (m > 0) {
14047c6ae99SBarry Smith       n = csize/m;
14147c6ae99SBarry Smith       if (m*n == csize) break;
14247c6ae99SBarry Smith       m--;
14347c6ae99SBarry Smith     }
14447c6ae99SBarry Smith     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
14547c6ae99SBarry Smith 
14647c6ae99SBarry Smith     x = M/m + ((M % m) > ((csize-1) % m));
14747c6ae99SBarry Smith     y = (N + (csize-1)/m)/n;
14847c6ae99SBarry Smith   } while ((x < 4 || y < 4) && csize > 1);
14947c6ae99SBarry Smith   if (size != csize) {
15047c6ae99SBarry Smith     MPI_Group    entire_group,sub_group;
15147c6ae99SBarry Smith     PetscMPIInt  i,*groupies;
15247c6ae99SBarry Smith 
15347c6ae99SBarry Smith     ierr = MPI_Comm_group(comm,&entire_group);CHKERRQ(ierr);
15447c6ae99SBarry Smith     ierr = PetscMalloc(csize*sizeof(PetscInt),&groupies);CHKERRQ(ierr);
15547c6ae99SBarry Smith     for (i=0; i<csize; i++) {
15647c6ae99SBarry Smith       groupies[i] = (rank/csize)*csize + i;
15747c6ae99SBarry Smith     }
15847c6ae99SBarry Smith     ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);CHKERRQ(ierr);
15947c6ae99SBarry Smith     ierr = PetscFree(groupies);CHKERRQ(ierr);
16047c6ae99SBarry Smith     ierr = MPI_Comm_create(comm,sub_group,outcomm);CHKERRQ(ierr);
16147c6ae99SBarry Smith     ierr = MPI_Group_free(&entire_group);CHKERRQ(ierr);
16247c6ae99SBarry Smith     ierr = MPI_Group_free(&sub_group);CHKERRQ(ierr);
163aa219208SBarry Smith     ierr = PetscInfo1(0,"DMDASplitComm2d:Creating redundant coarse problems of size %d\n",csize);CHKERRQ(ierr);
16447c6ae99SBarry Smith   } else {
16547c6ae99SBarry Smith     *outcomm = comm;
16647c6ae99SBarry Smith   }
16747c6ae99SBarry Smith   PetscFunctionReturn(0);
16847c6ae99SBarry Smith }
16947c6ae99SBarry Smith 
17047c6ae99SBarry Smith #if defined(new)
17147c6ae99SBarry Smith #undef __FUNCT__
172aa219208SBarry Smith #define __FUNCT__ "DMDAGetDiagonal_MFFD"
17347c6ae99SBarry Smith /*
174aa219208SBarry Smith   DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
175aa219208SBarry Smith     function lives on a DMDA
17647c6ae99SBarry Smith 
17747c6ae99SBarry Smith         y ~= (F(u + ha) - F(u))/h,
17847c6ae99SBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
17947c6ae99SBarry Smith         u = current iterate
18047c6ae99SBarry Smith         h = difference interval
18147c6ae99SBarry Smith */
182aa219208SBarry Smith PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a)
18347c6ae99SBarry Smith {
18447c6ae99SBarry Smith   PetscScalar    h,*aa,*ww,v;
18547c6ae99SBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
18647c6ae99SBarry Smith   PetscErrorCode ierr;
18747c6ae99SBarry Smith   PetscInt       gI,nI;
18847c6ae99SBarry Smith   MatStencil     stencil;
189aa219208SBarry Smith   DMDALocalInfo  info;
19047c6ae99SBarry Smith 
19147c6ae99SBarry Smith   PetscFunctionBegin;
19247c6ae99SBarry Smith   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
19347c6ae99SBarry Smith   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
19447c6ae99SBarry Smith 
19547c6ae99SBarry Smith   ierr = VecGetArray(U,&ww);CHKERRQ(ierr);
19647c6ae99SBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
19747c6ae99SBarry Smith 
19847c6ae99SBarry Smith   nI = 0;
19947c6ae99SBarry Smith     h  = ww[gI];
20047c6ae99SBarry Smith     if (h == 0.0) h = 1.0;
20147c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX)
20247c6ae99SBarry Smith     if (h < umin && h >= 0.0)      h = umin;
20347c6ae99SBarry Smith     else if (h < 0.0 && h > -umin) h = -umin;
20447c6ae99SBarry Smith #else
20547c6ae99SBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
20647c6ae99SBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
20747c6ae99SBarry Smith #endif
20847c6ae99SBarry Smith     h     *= epsilon;
20947c6ae99SBarry Smith 
21047c6ae99SBarry Smith     ww[gI] += h;
21147c6ae99SBarry Smith     ierr          = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
21247c6ae99SBarry Smith     aa[nI]  = (v - aa[nI])/h;
21347c6ae99SBarry Smith     ww[gI] -= h;
21447c6ae99SBarry Smith     nI++;
21547c6ae99SBarry Smith   }
21647c6ae99SBarry Smith   ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr);
21747c6ae99SBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
21847c6ae99SBarry Smith   PetscFunctionReturn(0);
21947c6ae99SBarry Smith }
22047c6ae99SBarry Smith #endif
22147c6ae99SBarry Smith 
22247c6ae99SBarry Smith #undef __FUNCT__
2239a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_2D"
2247087cfbeSBarry Smith PetscErrorCode  DMSetUp_DA_2D(DM da)
22547c6ae99SBarry Smith {
22647c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
22747c6ae99SBarry Smith   const PetscInt   M            = dd->M;
22847c6ae99SBarry Smith   const PetscInt   N            = dd->N;
22947c6ae99SBarry Smith   PetscInt         m            = dd->m;
23047c6ae99SBarry Smith   PetscInt         n            = dd->n;
23188661749SPeter Brune   PetscInt         o            = dd->overlap;
23247c6ae99SBarry Smith   const PetscInt   dof          = dd->w;
23347c6ae99SBarry Smith   const PetscInt   s            = dd->s;
23419fd82e9SBarry Smith   DMDABoundaryType bx         = dd->bx;
23519fd82e9SBarry Smith   DMDABoundaryType by         = dd->by;
23619fd82e9SBarry Smith   DMDAStencilType  stencil_type = dd->stencil_type;
23747c6ae99SBarry Smith   PetscInt         *lx           = dd->lx;
23847c6ae99SBarry Smith   PetscInt         *ly           = dd->ly;
23947c6ae99SBarry Smith   MPI_Comm         comm;
24047c6ae99SBarry Smith   PetscMPIInt      rank,size;
241ce00eea3SSatish Balay   PetscInt         xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,start,end,IXs,IXe,IYs,IYe;
242ce00eea3SSatish Balay   PetscInt         up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn,*idx_cpy;
243ce00eea3SSatish Balay   const PetscInt   *idx_full;
244db87c5ecSEthan Coon   PetscInt         xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
24547c6ae99SBarry Smith   PetscInt         s_x,s_y; /* s proportionalized to w */
24647c6ae99SBarry Smith   PetscInt         sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
24747c6ae99SBarry Smith   Vec              local,global;
24847c6ae99SBarry Smith   VecScatter       ltog,gtol;
249ce00eea3SSatish Balay   IS               to,from,ltogis;
25047c6ae99SBarry Smith   PetscErrorCode   ierr;
25147c6ae99SBarry Smith 
25247c6ae99SBarry Smith   PetscFunctionBegin;
2533855c12bSBarry Smith   if (dof < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
2543855c12bSBarry Smith   if (s < 0) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
25547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
2563855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES)
2573855c12bSBarry Smith   if (((Petsc64bitInt) M)*((Petsc64bitInt) N)*((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);
2583855c12bSBarry Smith #endif
25947c6ae99SBarry Smith 
26047c6ae99SBarry Smith   if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
26147c6ae99SBarry Smith   if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
26247c6ae99SBarry Smith 
26347c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
26447c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
26547c6ae99SBarry Smith 
26647c6ae99SBarry Smith   dd->dim         = 2;
26747c6ae99SBarry Smith   ierr = PetscMalloc(dof*sizeof(char*),&dd->fieldname);CHKERRQ(ierr);
26847c6ae99SBarry Smith   ierr = PetscMemzero(dd->fieldname,dof*sizeof(char*));CHKERRQ(ierr);
26947c6ae99SBarry Smith 
27047c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
27147c6ae99SBarry Smith     if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
27247c6ae99SBarry Smith     else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
27347c6ae99SBarry Smith   }
27447c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
27547c6ae99SBarry Smith     if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
27647c6ae99SBarry Smith     else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
27747c6ae99SBarry Smith   }
27847c6ae99SBarry Smith 
27947c6ae99SBarry Smith   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
28047c6ae99SBarry Smith     if (n != PETSC_DECIDE) {
28147c6ae99SBarry Smith       m = size/n;
28247c6ae99SBarry Smith     } else if (m != PETSC_DECIDE) {
28347c6ae99SBarry Smith       n = size/m;
28447c6ae99SBarry Smith     } else {
28547c6ae99SBarry Smith       /* try for squarish distribution */
28647c6ae99SBarry Smith       m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)N)));
28747c6ae99SBarry Smith       if (!m) m = 1;
28847c6ae99SBarry Smith       while (m > 0) {
28947c6ae99SBarry Smith 	n = size/m;
29047c6ae99SBarry Smith 	if (m*n == size) break;
29147c6ae99SBarry Smith 	m--;
29247c6ae99SBarry Smith       }
29347c6ae99SBarry Smith       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
29447c6ae99SBarry Smith     }
29547c6ae99SBarry Smith     if (m*n != size) SETERRQ(comm,PETSC_ERR_PLIB,"Unable to create partition, check the size of the communicator and input m and n ");
29647c6ae99SBarry Smith   } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
29747c6ae99SBarry Smith 
29847c6ae99SBarry Smith   if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
29947c6ae99SBarry Smith   if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
30047c6ae99SBarry Smith 
30147c6ae99SBarry Smith   /*
30247c6ae99SBarry Smith      Determine locally owned region
30347c6ae99SBarry Smith      xs is the first local node number, x is the number of local nodes
30447c6ae99SBarry Smith   */
30547c6ae99SBarry Smith   if (!lx) {
30647c6ae99SBarry Smith     ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
30747c6ae99SBarry Smith     lx = dd->lx;
30847c6ae99SBarry Smith     for (i=0; i<m; i++) {
30947c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > i);
31047c6ae99SBarry Smith     }
31147c6ae99SBarry Smith   }
31247c6ae99SBarry Smith   x  = lx[rank % m];
31347c6ae99SBarry Smith   xs = 0;
31447c6ae99SBarry Smith   for (i=0; i<(rank % m); i++) {
31547c6ae99SBarry Smith     xs += lx[i];
31647c6ae99SBarry Smith   }
31747c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
31847c6ae99SBarry Smith   left = xs;
31947c6ae99SBarry Smith   for (i=(rank % m); i<m; i++) {
32047c6ae99SBarry Smith     left += lx[i];
32147c6ae99SBarry Smith   }
32247c6ae99SBarry Smith   if (left != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M: %D %D",left,M);
32347c6ae99SBarry Smith #endif
32447c6ae99SBarry Smith 
32547c6ae99SBarry Smith   /*
32647c6ae99SBarry Smith      Determine locally owned region
32747c6ae99SBarry Smith      ys is the first local node number, y is the number of local nodes
32847c6ae99SBarry Smith   */
32947c6ae99SBarry Smith   if (!ly) {
33047c6ae99SBarry Smith     ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
33147c6ae99SBarry Smith     ly = dd->ly;
33247c6ae99SBarry Smith     for (i=0; i<n; i++) {
33347c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > i);
33447c6ae99SBarry Smith     }
33547c6ae99SBarry Smith   }
33647c6ae99SBarry Smith   y  = ly[rank/m];
33747c6ae99SBarry Smith   ys = 0;
33847c6ae99SBarry Smith   for (i=0; i<(rank/m); i++) {
33947c6ae99SBarry Smith     ys += ly[i];
34047c6ae99SBarry Smith   }
34147c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
34247c6ae99SBarry Smith   left = ys;
34347c6ae99SBarry Smith   for (i=(rank/m); i<n; i++) {
34447c6ae99SBarry Smith     left += ly[i];
34547c6ae99SBarry Smith   }
34647c6ae99SBarry Smith   if (left != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of ly across processors not equal to N: %D %D",left,N);
34747c6ae99SBarry Smith #endif
34847c6ae99SBarry Smith 
349bcea557cSEthan Coon   /*
350bcea557cSEthan Coon    check if the scatter requires more than one process neighbor or wraps around
351bcea557cSEthan Coon    the domain more than once
352bcea557cSEthan Coon   */
35388661749SPeter Brune   if ((x < s+o) && ((m > 1) || (bx == DMDA_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s+o);
35488661749SPeter Brune   if ((y < s+o) && ((n > 1) || (by == DMDA_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s+o);
35547c6ae99SBarry Smith   xe = xs + x;
35647c6ae99SBarry Smith   ye = ys + y;
35747c6ae99SBarry Smith 
358ce00eea3SSatish Balay   /* determine ghost region (Xs) and region scattered into (IXs)  */
35988661749SPeter Brune   if (xs-s-o > 0) {
36088661749SPeter Brune     Xs = xs - s - o; IXs = xs - s - o;
36188661749SPeter Brune   } else {
36288661749SPeter Brune     if (bx) {
36388661749SPeter Brune       Xs = xs - s;
36488661749SPeter Brune     } else {
36588661749SPeter Brune       Xs = 0;
36688661749SPeter Brune     }
36788661749SPeter Brune     IXs = 0;
36888661749SPeter Brune   }
36988661749SPeter Brune   if (xe+s+o <= M) {
37088661749SPeter Brune     Xe = xe + s + o; IXe = xe + s + o;
37188661749SPeter Brune   } else {
37288661749SPeter Brune     if (bx) {
37388661749SPeter Brune       Xs = xs - s - o; Xe = xe + s;
37488661749SPeter Brune     } else {
37588661749SPeter Brune       Xe = M;
37688661749SPeter Brune     }
37788661749SPeter Brune     IXe = M;
37888661749SPeter Brune   }
37947c6ae99SBarry Smith 
38088661749SPeter Brune   if (bx == DMDA_BOUNDARY_PERIODIC) {
38188661749SPeter Brune     IXs = xs - s - o;
38288661749SPeter Brune     IXe = xe + s + o;
38388661749SPeter Brune     Xs = xs - s - o;
38488661749SPeter Brune     Xe = xe + s + o;
38588661749SPeter Brune   }
38647c6ae99SBarry Smith 
38788661749SPeter Brune   if (ys-s-o > 0) {
38888661749SPeter Brune     Ys = ys - s - o; IYs = ys - s - o;
38988661749SPeter Brune   } else {
39088661749SPeter Brune     if (by) {
39188661749SPeter Brune       Ys = ys - s;
39288661749SPeter Brune     } else {
39388661749SPeter Brune       Ys = 0;
39488661749SPeter Brune     }
39588661749SPeter Brune     IYs = 0;
39688661749SPeter Brune   }
39788661749SPeter Brune   if (ye+s+o <= N) {
39888661749SPeter Brune     Ye = ye + s + o; IYe = ye + s + o;
39988661749SPeter Brune   } else {
40088661749SPeter Brune     if (by) {
40188661749SPeter Brune       Ye = ye + s;
40288661749SPeter Brune     } else {
40388661749SPeter Brune       Ye = N;
40488661749SPeter Brune     }
40588661749SPeter Brune     IYe = N;
40688661749SPeter Brune   }
40788661749SPeter Brune 
40888661749SPeter Brune   if (by == DMDA_BOUNDARY_PERIODIC) {
40988661749SPeter Brune     IYs = ys - s - o;
41088661749SPeter Brune     IYe = ye + s + o;
41188661749SPeter Brune     Ys = ys - s - o;
41288661749SPeter Brune     Ye = ye + s + o;
41388661749SPeter Brune   }
41488661749SPeter Brune 
41588661749SPeter Brune   /* stencil length in each direction */
41688661749SPeter Brune   s_x = s + o;
41788661749SPeter Brune   s_y = s + o;
41847c6ae99SBarry Smith 
41947c6ae99SBarry Smith   /* determine starting point of each processor */
42047c6ae99SBarry Smith   nn    = x*y;
42147c6ae99SBarry Smith   ierr  = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr);
42247c6ae99SBarry Smith   ierr  = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
42347c6ae99SBarry Smith   bases[0] = 0;
42447c6ae99SBarry Smith   for (i=1; i<=size; i++) {
42547c6ae99SBarry Smith     bases[i] = ldims[i-1];
42647c6ae99SBarry Smith   }
42747c6ae99SBarry Smith   for (i=1; i<=size; i++) {
42847c6ae99SBarry Smith     bases[i] += bases[i-1];
42947c6ae99SBarry Smith   }
430ce00eea3SSatish Balay   base = bases[rank]*dof;
43147c6ae99SBarry Smith 
43247c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
433ce00eea3SSatish Balay   dd->Nlocal = x*y*dof;
434778a2246SBarry Smith   ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr);
435ce00eea3SSatish Balay   dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof;
436778a2246SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,0,&local);CHKERRQ(ierr);
43747c6ae99SBarry Smith 
43847c6ae99SBarry Smith   /* generate appropriate vector scatters */
43947c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
44047c6ae99SBarry Smith   ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr);
441ce00eea3SSatish Balay   ierr = ISCreateStride(comm,x*y*dof,start,1,&to);CHKERRQ(ierr);
44247c6ae99SBarry Smith 
443db87c5ecSEthan Coon   count = x*y;
444ce00eea3SSatish Balay   ierr = PetscMalloc(x*y*sizeof(PetscInt),&idx);CHKERRQ(ierr);
445ce00eea3SSatish Balay   left = xs - Xs; right = left + x;
446ce00eea3SSatish Balay   down = ys - Ys; up = down + y;
44747c6ae99SBarry Smith   count = 0;
44847c6ae99SBarry Smith   for (i=down; i<up; i++) {
449ce00eea3SSatish Balay     for (j=left; j<right; j++) {
450ce00eea3SSatish Balay       idx[count++] = i*(Xe-Xs) + j;
45147c6ae99SBarry Smith     }
45247c6ae99SBarry Smith   }
45347c6ae99SBarry Smith 
454ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
45547c6ae99SBarry Smith   ierr = VecScatterCreate(local,from,global,to,&ltog);CHKERRQ(ierr);
45647c6ae99SBarry Smith   ierr = PetscLogObjectParent(dd,ltog);CHKERRQ(ierr);
457fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
458fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
45947c6ae99SBarry Smith 
460ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
461ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
46288661749SPeter Brune   if (stencil_type == DMDA_STENCIL_BOX || o > 0) {
463db87c5ecSEthan Coon     count = (IXe-IXs)*(IYe-IYs);
464db87c5ecSEthan Coon     ierr  = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
465ce00eea3SSatish Balay 
466ce00eea3SSatish Balay     left = IXs - Xs; right = left + (IXe-IXs);
467ce00eea3SSatish Balay     down = IYs - Ys; up = down + (IYe-IYs);
468ce00eea3SSatish Balay     count = 0;
469ce00eea3SSatish Balay     for (i=down; i<up; i++) {
470ce00eea3SSatish Balay       for (j=left; j<right; j++) {
471ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
472ce00eea3SSatish Balay       }
473ce00eea3SSatish Balay     }
474ce00eea3SSatish Balay     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
475ce00eea3SSatish Balay 
47647c6ae99SBarry Smith   } else {
47747c6ae99SBarry Smith     /* must drop into cross shape region */
47847c6ae99SBarry Smith     /*       ---------|
47947c6ae99SBarry Smith             |  top    |
480ce00eea3SSatish Balay          |---         ---| up
48147c6ae99SBarry Smith          |   middle      |
48247c6ae99SBarry Smith          |               |
483ce00eea3SSatish Balay          ----         ---- down
48447c6ae99SBarry Smith             | bottom  |
48547c6ae99SBarry Smith             -----------
48647c6ae99SBarry Smith          Xs xs        xe Xe */
487db87c5ecSEthan Coon     count = (ys-IYs)*x + y*(IXe-IXs) + (IYe-ye)*x;
488db87c5ecSEthan Coon     ierr  = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
489ce00eea3SSatish Balay 
490ce00eea3SSatish Balay     left = xs - Xs; right = left + x;
491ce00eea3SSatish Balay     down = ys - Ys; up = down + y;
49247c6ae99SBarry Smith     count = 0;
493ce00eea3SSatish Balay     /* bottom */
494ce00eea3SSatish Balay     for (i=(IYs-Ys); i<down; i++) {
495ce00eea3SSatish Balay       for (j=left; j<right; j++) {
496ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
49747c6ae99SBarry Smith       }
49847c6ae99SBarry Smith     }
49947c6ae99SBarry Smith     /* middle */
50047c6ae99SBarry Smith     for (i=down; i<up; i++) {
501ce00eea3SSatish Balay       for (j=(IXs-Xs); j<(IXe-Xs); j++) {
502ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
50347c6ae99SBarry Smith       }
50447c6ae99SBarry Smith     }
50547c6ae99SBarry Smith     /* top */
506ce00eea3SSatish Balay     for (i=up; i<up+IYe-ye; i++) {
507ce00eea3SSatish Balay       for (j=left; j<right; j++) {
508ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
50947c6ae99SBarry Smith       }
51047c6ae99SBarry Smith     }
51147c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
51247c6ae99SBarry Smith   }
51347c6ae99SBarry Smith 
51447c6ae99SBarry Smith 
51547c6ae99SBarry Smith   /* determine who lies on each side of us stored in    n6 n7 n8
51647c6ae99SBarry Smith                                                         n3    n5
51747c6ae99SBarry Smith                                                         n0 n1 n2
51847c6ae99SBarry Smith   */
51947c6ae99SBarry Smith 
52047c6ae99SBarry Smith   /* Assume the Non-Periodic Case */
52147c6ae99SBarry Smith   n1 = rank - m;
52247c6ae99SBarry Smith   if (rank % m) {
52347c6ae99SBarry Smith     n0 = n1 - 1;
52447c6ae99SBarry Smith   } else {
52547c6ae99SBarry Smith     n0 = -1;
52647c6ae99SBarry Smith   }
52747c6ae99SBarry Smith   if ((rank+1) % m) {
52847c6ae99SBarry Smith     n2 = n1 + 1;
52947c6ae99SBarry Smith     n5 = rank + 1;
53047c6ae99SBarry Smith     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
53147c6ae99SBarry Smith   } else {
53247c6ae99SBarry Smith     n2 = -1; n5 = -1; n8 = -1;
53347c6ae99SBarry Smith   }
53447c6ae99SBarry Smith   if (rank % m) {
53547c6ae99SBarry Smith     n3 = rank - 1;
53647c6ae99SBarry Smith     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
53747c6ae99SBarry Smith   } else {
53847c6ae99SBarry Smith     n3 = -1; n6 = -1;
53947c6ae99SBarry Smith   }
54047c6ae99SBarry Smith   n7 = rank + m; if (n7 >= m*n) n7 = -1;
54147c6ae99SBarry Smith 
5421321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC && by == DMDA_BOUNDARY_PERIODIC) {
54347c6ae99SBarry Smith   /* Modify for Periodic Cases */
54447c6ae99SBarry Smith     /* Handle all four corners */
54547c6ae99SBarry Smith     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
54647c6ae99SBarry Smith     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
54747c6ae99SBarry Smith     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
54847c6ae99SBarry Smith     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
54947c6ae99SBarry Smith 
55047c6ae99SBarry Smith     /* Handle Top and Bottom Sides */
55147c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
55247c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
55347c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
55447c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
55547c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
55647c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
55747c6ae99SBarry Smith 
55847c6ae99SBarry Smith     /* Handle Left and Right Sides */
55947c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
56047c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
56147c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
56247c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
56347c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
56447c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
5651321219cSEthan Coon   } else if (by == DMDA_BOUNDARY_PERIODIC) {  /* Handle Top and Bottom Sides */
566ce00eea3SSatish Balay     if (n1 < 0) n1 = rank + m * (n-1);
567ce00eea3SSatish Balay     if (n7 < 0) n7 = rank - m * (n-1);
568ce00eea3SSatish Balay     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
569ce00eea3SSatish Balay     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
570ce00eea3SSatish Balay     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
571ce00eea3SSatish Balay     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
5721321219cSEthan Coon   } else if (bx == DMDA_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */
573ce00eea3SSatish Balay     if (n3 < 0) n3 = rank + (m-1);
574ce00eea3SSatish Balay     if (n5 < 0) n5 = rank - (m-1);
575ce00eea3SSatish Balay     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
576ce00eea3SSatish Balay     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
577ce00eea3SSatish Balay     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
578ce00eea3SSatish Balay     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
57947c6ae99SBarry Smith   }
580ce00eea3SSatish Balay 
58147c6ae99SBarry Smith   ierr = PetscMalloc(9*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr);
58247c6ae99SBarry Smith   dd->neighbors[0] = n0;
58347c6ae99SBarry Smith   dd->neighbors[1] = n1;
58447c6ae99SBarry Smith   dd->neighbors[2] = n2;
58547c6ae99SBarry Smith   dd->neighbors[3] = n3;
58647c6ae99SBarry Smith   dd->neighbors[4] = rank;
58747c6ae99SBarry Smith   dd->neighbors[5] = n5;
58847c6ae99SBarry Smith   dd->neighbors[6] = n6;
58947c6ae99SBarry Smith   dd->neighbors[7] = n7;
59047c6ae99SBarry Smith   dd->neighbors[8] = n8;
59147c6ae99SBarry Smith 
59288661749SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR && o == 0) {
59347c6ae99SBarry Smith     /* save corner processor numbers */
59447c6ae99SBarry Smith     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
59547c6ae99SBarry Smith     n0 = n2 = n6 = n8 = -1;
59647c6ae99SBarry Smith   }
59747c6ae99SBarry Smith 
598ce00eea3SSatish Balay   ierr = PetscMalloc((Xe-Xs)*(Ye-Ys)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
599ce00eea3SSatish Balay   ierr = PetscLogObjectMemory(da,(Xe-Xs)*(Ye-Ys)*sizeof(PetscInt));CHKERRQ(ierr);
60047c6ae99SBarry Smith 
601ce00eea3SSatish Balay   nn = 0;
60247c6ae99SBarry Smith   xbase = bases[rank];
60347c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
60447c6ae99SBarry Smith     if (n0 >= 0) { /* left below */
605ce00eea3SSatish Balay       x_t = lx[n0 % m];
60647c6ae99SBarry Smith       y_t = ly[(n0/m)];
60747c6ae99SBarry Smith       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
60847c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
60947c6ae99SBarry Smith     }
61047c6ae99SBarry Smith     if (n1 >= 0) { /* directly below */
61147c6ae99SBarry Smith       x_t = x;
61247c6ae99SBarry Smith       y_t = ly[(n1/m)];
61347c6ae99SBarry Smith       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
61447c6ae99SBarry Smith       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
61547c6ae99SBarry Smith     }
61647c6ae99SBarry Smith     if (n2 >= 0) { /* right below */
617ce00eea3SSatish Balay       x_t = lx[n2 % m];
61847c6ae99SBarry Smith       y_t = ly[(n2/m)];
61947c6ae99SBarry Smith       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
62047c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
62147c6ae99SBarry Smith     }
62247c6ae99SBarry Smith   }
62347c6ae99SBarry Smith 
62447c6ae99SBarry Smith   for (i=0; i<y; i++) {
62547c6ae99SBarry Smith     if (n3 >= 0) { /* directly left */
626ce00eea3SSatish Balay       x_t = lx[n3 % m];
62747c6ae99SBarry Smith       /* y_t = y; */
62847c6ae99SBarry Smith       s_t = bases[n3] + (i+1)*x_t - s_x;
62947c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
63047c6ae99SBarry Smith     }
63147c6ae99SBarry Smith 
63247c6ae99SBarry Smith     for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
63347c6ae99SBarry Smith 
63447c6ae99SBarry Smith     if (n5 >= 0) { /* directly right */
635ce00eea3SSatish Balay       x_t = lx[n5 % m];
63647c6ae99SBarry Smith       /* y_t = y; */
63747c6ae99SBarry Smith       s_t = bases[n5] + (i)*x_t;
63847c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
63947c6ae99SBarry Smith     }
64047c6ae99SBarry Smith   }
64147c6ae99SBarry Smith 
64247c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
64347c6ae99SBarry Smith     if (n6 >= 0) { /* left above */
644ce00eea3SSatish Balay       x_t = lx[n6 % m];
64547c6ae99SBarry Smith       /* y_t = ly[(n6/m)]; */
64647c6ae99SBarry Smith       s_t = bases[n6] + (i)*x_t - s_x;
64747c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
64847c6ae99SBarry Smith     }
64947c6ae99SBarry Smith     if (n7 >= 0) { /* directly above */
65047c6ae99SBarry Smith       x_t = x;
65147c6ae99SBarry Smith       /* y_t = ly[(n7/m)]; */
65247c6ae99SBarry Smith       s_t = bases[n7] + (i-1)*x_t;
65347c6ae99SBarry Smith       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
65447c6ae99SBarry Smith     }
65547c6ae99SBarry Smith     if (n8 >= 0) { /* right above */
656ce00eea3SSatish Balay       x_t = lx[n8 % m];
65747c6ae99SBarry Smith       /* y_t = ly[(n8/m)]; */
65847c6ae99SBarry Smith       s_t = bases[n8] + (i-1)*x_t;
65947c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
66047c6ae99SBarry Smith     }
66147c6ae99SBarry Smith   }
66247c6ae99SBarry Smith 
663ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
66447c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
66547c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr);
666fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
667fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
66847c6ae99SBarry Smith 
66988661749SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR && o == 0) {
670ce00eea3SSatish Balay     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
671ce00eea3SSatish Balay   }
672ce00eea3SSatish Balay 
67388661749SPeter Brune   if (((stencil_type == DMDA_STENCIL_STAR) ||
6741321219cSEthan Coon        (bx && bx != DMDA_BOUNDARY_PERIODIC) ||
67588661749SPeter Brune        (by && by != DMDA_BOUNDARY_PERIODIC)) && o == 0) {
67647c6ae99SBarry Smith     /*
67747c6ae99SBarry Smith         Recompute the local to global mappings, this time keeping the
678ce00eea3SSatish Balay       information about the cross corner processor numbers and any ghosted
679ce00eea3SSatish Balay       but not periodic indices.
68047c6ae99SBarry Smith     */
68147c6ae99SBarry Smith     nn = 0;
68247c6ae99SBarry Smith     xbase = bases[rank];
68347c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
68447c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
685ce00eea3SSatish Balay         x_t = lx[n0 % m];
68647c6ae99SBarry Smith         y_t = ly[(n0/m)];
68747c6ae99SBarry Smith         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
68847c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
689ce00eea3SSatish Balay       } else if (xs-Xs > 0 && ys-Ys > 0) {
690ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
69147c6ae99SBarry Smith       }
69247c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
69347c6ae99SBarry Smith         x_t = x;
69447c6ae99SBarry Smith         y_t = ly[(n1/m)];
69547c6ae99SBarry Smith         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
69647c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
697ce00eea3SSatish Balay       } else if (ys-Ys > 0) {
698ce00eea3SSatish Balay         for (j=0; j<x; j++) { idx[nn++] = -1;}
69947c6ae99SBarry Smith       }
70047c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
701ce00eea3SSatish Balay         x_t = lx[n2 % m];
70247c6ae99SBarry Smith         y_t = ly[(n2/m)];
70347c6ae99SBarry Smith         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
70447c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
705ce00eea3SSatish Balay       } else if (Xe-xe> 0 && ys-Ys > 0) {
706ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
70747c6ae99SBarry Smith       }
70847c6ae99SBarry Smith     }
70947c6ae99SBarry Smith 
71047c6ae99SBarry Smith     for (i=0; i<y; i++) {
71147c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
712ce00eea3SSatish Balay         x_t = lx[n3 % m];
71347c6ae99SBarry Smith         /* y_t = y; */
71447c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x;
71547c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
716ce00eea3SSatish Balay       } else if (xs-Xs > 0) {
717ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
71847c6ae99SBarry Smith       }
71947c6ae99SBarry Smith 
72047c6ae99SBarry Smith       for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
72147c6ae99SBarry Smith 
72247c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
723ce00eea3SSatish Balay         x_t = lx[n5 % m];
72447c6ae99SBarry Smith         /* y_t = y; */
72547c6ae99SBarry Smith         s_t = bases[n5] + (i)*x_t;
72647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
727ce00eea3SSatish Balay       } else if (Xe-xe > 0) {
728ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
72947c6ae99SBarry Smith       }
73047c6ae99SBarry Smith     }
73147c6ae99SBarry Smith 
73247c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
73347c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
734ce00eea3SSatish Balay         x_t = lx[n6 % m];
73547c6ae99SBarry Smith         /* y_t = ly[(n6/m)]; */
73647c6ae99SBarry Smith         s_t = bases[n6] + (i)*x_t - s_x;
73747c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
738ce00eea3SSatish Balay       } else if (xs-Xs > 0 && Ye-ye > 0) {
739ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
74047c6ae99SBarry Smith       }
74147c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
74247c6ae99SBarry Smith         x_t = x;
74347c6ae99SBarry Smith         /* y_t = ly[(n7/m)]; */
74447c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t;
74547c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
746ce00eea3SSatish Balay       } else if (Ye-ye > 0) {
747ce00eea3SSatish Balay         for (j=0; j<x; j++) { idx[nn++] = -1;}
74847c6ae99SBarry Smith       }
74947c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
750ce00eea3SSatish Balay         x_t = lx[n8 % m];
75147c6ae99SBarry Smith         /* y_t = ly[(n8/m)]; */
75247c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t;
75347c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
754ce00eea3SSatish Balay       } else if (Xe-xe > 0 && Ye-ye > 0) {
755ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
75647c6ae99SBarry Smith       }
75747c6ae99SBarry Smith     }
75847c6ae99SBarry Smith   }
759ce00eea3SSatish Balay   /*
760ce00eea3SSatish Balay      Set the local to global ordering in the global vector, this allows use
761ce00eea3SSatish Balay      of VecSetValuesLocal().
762ce00eea3SSatish Balay   */
763ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_OWN_POINTER,&ltogis);CHKERRQ(ierr);
764ce00eea3SSatish Balay   ierr = PetscMalloc(nn*dof*sizeof(PetscInt),&idx_cpy);CHKERRQ(ierr);
765db87c5ecSEthan Coon   ierr = PetscLogObjectMemory(da,nn*dof*sizeof(PetscInt));CHKERRQ(ierr);
766ce00eea3SSatish Balay   ierr = ISGetIndices(ltogis, &idx_full);
767ce00eea3SSatish Balay   ierr = PetscMemcpy(idx_cpy,idx_full,nn*dof*sizeof(PetscInt));CHKERRQ(ierr);
768ce00eea3SSatish Balay   ierr = ISRestoreIndices(ltogis, &idx_full);
769ce00eea3SSatish Balay   ierr = ISLocalToGlobalMappingCreateIS(ltogis,&da->ltogmap);CHKERRQ(ierr);
770ce00eea3SSatish Balay   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
771fcfd50ebSBarry Smith   ierr = ISDestroy(&ltogis);CHKERRQ(ierr);
772ce00eea3SSatish Balay   ierr = ISLocalToGlobalMappingBlock(da->ltogmap,dd->w,&da->ltogmapb);CHKERRQ(ierr);
773ce00eea3SSatish Balay   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
77447c6ae99SBarry Smith 
775ce00eea3SSatish Balay   ierr = PetscFree2(bases,ldims);CHKERRQ(ierr);
77647c6ae99SBarry Smith   dd->m  = m;  dd->n  = n;
777ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
778ce00eea3SSatish Balay   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
779ce00eea3SSatish Balay   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
78047c6ae99SBarry Smith 
781fcfd50ebSBarry Smith   ierr = VecDestroy(&local);CHKERRQ(ierr);
782fcfd50ebSBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
78347c6ae99SBarry Smith 
78447c6ae99SBarry Smith   dd->gtol      = gtol;
78547c6ae99SBarry Smith   dd->ltog      = ltog;
786ce00eea3SSatish Balay   dd->idx       = idx_cpy;
787ce00eea3SSatish Balay   dd->Nl        = nn*dof;
78847c6ae99SBarry Smith   dd->base      = base;
7899a42bb27SBarry Smith   da->ops->view = DMView_DA_2d;
79047c6ae99SBarry Smith   dd->ltol = PETSC_NULL;
79147c6ae99SBarry Smith   dd->ao   = PETSC_NULL;
79247c6ae99SBarry Smith 
79347c6ae99SBarry Smith   PetscFunctionReturn(0);
79447c6ae99SBarry Smith }
79547c6ae99SBarry Smith 
79647c6ae99SBarry Smith #undef __FUNCT__
797aa219208SBarry Smith #define __FUNCT__ "DMDACreate2d"
79847c6ae99SBarry Smith /*@C
799aa219208SBarry Smith    DMDACreate2d -  Creates an object that will manage the communication of  two-dimensional
80047c6ae99SBarry Smith    regular array data that is distributed across some processors.
80147c6ae99SBarry Smith 
80247c6ae99SBarry Smith    Collective on MPI_Comm
80347c6ae99SBarry Smith 
80447c6ae99SBarry Smith    Input Parameters:
80547c6ae99SBarry Smith +  comm - MPI communicator
8061321219cSEthan Coon .  bx,by - type of ghost nodes the array have.
8071321219cSEthan Coon          Use one of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC.
808aa219208SBarry Smith .  stencil_type - stencil type.  Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
80947c6ae99SBarry Smith .  M,N - global dimension in each direction of the array (use -M and or -N to indicate that it may be set to a different value
81047c6ae99SBarry Smith             from the command line with -da_grid_x <M> -da_grid_y <N>)
81147c6ae99SBarry Smith .  m,n - corresponding number of processors in each dimension
81247c6ae99SBarry Smith          (or PETSC_DECIDE to have calculated)
81347c6ae99SBarry Smith .  dof - number of degrees of freedom per node
81447c6ae99SBarry Smith .  s - stencil width
81547c6ae99SBarry Smith -  lx, ly - arrays containing the number of nodes in each cell along
81647c6ae99SBarry Smith            the x and y coordinates, or PETSC_NULL. If non-null, these
81747c6ae99SBarry Smith            must be of length as m and n, and the corresponding
81847c6ae99SBarry Smith            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
81947c6ae99SBarry Smith            must be M, and the sum of the ly[] entries must be N.
82047c6ae99SBarry Smith 
82147c6ae99SBarry Smith    Output Parameter:
82247c6ae99SBarry Smith .  da - the resulting distributed array object
82347c6ae99SBarry Smith 
82447c6ae99SBarry Smith    Options Database Key:
825aa219208SBarry Smith +  -da_view - Calls DMView() at the conclusion of DMDACreate2d()
82647c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
82747c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
82847c6ae99SBarry Smith .  -da_processors_x <nx> - number of processors in x direction
82947c6ae99SBarry Smith .  -da_processors_y <ny> - number of processors in y direction
830e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement ratio in x direction
831e0f5d30fSBarry Smith .  -da_refine_y <ry> - refinement ratio in y direction
832e0f5d30fSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating, if M or N < 0
833e0f5d30fSBarry Smith 
83447c6ae99SBarry Smith 
83547c6ae99SBarry Smith    Level: beginner
83647c6ae99SBarry Smith 
83747c6ae99SBarry Smith    Notes:
838aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
839aa219208SBarry Smith    standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
84047c6ae99SBarry Smith    the standard 9-pt stencil.
84147c6ae99SBarry Smith 
842aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
843564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
844564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
84547c6ae99SBarry Smith 
84647c6ae99SBarry Smith .keywords: distributed array, create, two-dimensional
84747c6ae99SBarry Smith 
848aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
849aa219208SBarry Smith           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDASetRefinementFactor(),
850d461ba97SBarry Smith           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
85147c6ae99SBarry Smith 
85247c6ae99SBarry Smith @*/
853fe16a2e9SBarry Smith 
8541321219cSEthan Coon PetscErrorCode  DMDACreate2d(MPI_Comm comm,DMDABoundaryType bx,DMDABoundaryType by,DMDAStencilType stencil_type,
8559a42bb27SBarry Smith                           PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
85647c6ae99SBarry Smith {
85747c6ae99SBarry Smith   PetscErrorCode ierr;
85847c6ae99SBarry Smith 
85947c6ae99SBarry Smith   PetscFunctionBegin;
860aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
861aa219208SBarry Smith   ierr = DMDASetDim(*da, 2);CHKERRQ(ierr);
862aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, 1);CHKERRQ(ierr);
863aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr);
864755f258dSLisandro Dalcin   ierr = DMDASetBoundaryType(*da, bx, by, DMDA_BOUNDARY_NONE);CHKERRQ(ierr);
865aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
866aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
867aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
868aa219208SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, PETSC_NULL);CHKERRQ(ierr);
86947c6ae99SBarry Smith   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
8709a42bb27SBarry Smith   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
8719a42bb27SBarry Smith   ierr = DMSetUp(*da);CHKERRQ(ierr);
872*ca266f36SBarry Smith   ierr = DMViewFromOptions(*da,"-dm_view");CHKERRQ(ierr);
87347c6ae99SBarry Smith   PetscFunctionReturn(0);
87447c6ae99SBarry Smith }
875