xref: /petsc/src/dm/impls/da/da2.c (revision 832b7ceb95bd65cfd99922d89fc8cdbc69c04de1)
19a42bb27SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
39804daf3SBarry Smith #include <petscdraw.h>
447c6ae99SBarry Smith 
547c6ae99SBarry Smith #undef __FUNCT__
69a42bb27SBarry Smith #define __FUNCT__ "DMView_DA_2d"
79a42bb27SBarry Smith PetscErrorCode DMView_DA_2d(DM da,PetscViewer viewer)
847c6ae99SBarry Smith {
947c6ae99SBarry Smith   PetscErrorCode ierr;
1047c6ae99SBarry Smith   PetscMPIInt    rank;
119a42bb27SBarry Smith   PetscBool      iascii,isdraw,isbinary;
1247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
139a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
149a42bb27SBarry Smith   PetscBool ismatlab;
159a42bb27SBarry Smith #endif
1647c6ae99SBarry Smith 
1747c6ae99SBarry Smith   PetscFunctionBegin;
18ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr);
1947c6ae99SBarry Smith 
20251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
21251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
22251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
239a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
24251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr);
259a42bb27SBarry Smith #endif
2647c6ae99SBarry Smith   if (iascii) {
2747c6ae99SBarry Smith     PetscViewerFormat format;
2847c6ae99SBarry Smith 
2947c6ae99SBarry Smith     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
3047c6ae99SBarry Smith     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) {
31aa219208SBarry Smith       DMDALocalInfo info;
32aa219208SBarry Smith       ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
331575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
3447c6ae99SBarry 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);
3547c6ae99SBarry 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);
3647c6ae99SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
371575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
383da9ae13SJed Brown     } else {
393da9ae13SJed Brown       ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr);
4047c6ae99SBarry Smith     }
4147c6ae99SBarry Smith   } else if (isdraw) {
4247c6ae99SBarry Smith     PetscDraw      draw;
4347c6ae99SBarry Smith     double         ymin = -1*dd->s-1,ymax = dd->N+dd->s;
4447c6ae99SBarry Smith     double         xmin = -1*dd->s-1,xmax = dd->M+dd->s;
4547c6ae99SBarry Smith     double         x,y;
468ea3bf28SBarry Smith     PetscInt       base;
478ea3bf28SBarry Smith     const PetscInt *idx;
4847c6ae99SBarry Smith     char           node[10];
4947c6ae99SBarry Smith     PetscBool      isnull;
5047c6ae99SBarry Smith 
5147c6ae99SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
525b399a63SLisandro Dalcin     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
535b399a63SLisandro Dalcin     if (isnull) PetscFunctionReturn(0);
5447c6ae99SBarry Smith 
555b399a63SLisandro Dalcin     ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
565b399a63SLisandro Dalcin     ierr = PetscDrawClear(draw);CHKERRQ(ierr);
575b399a63SLisandro Dalcin     ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
585b399a63SLisandro Dalcin 
595b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
6047c6ae99SBarry Smith     /* first processor draw all node lines */
6147c6ae99SBarry Smith     if (!rank) {
6247c6ae99SBarry Smith       ymin = 0.0; ymax = dd->N - 1;
6347c6ae99SBarry Smith       for (xmin=0; xmin<dd->M; xmin++) {
6447c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
6547c6ae99SBarry Smith       }
6647c6ae99SBarry Smith       xmin = 0.0; xmax = dd->M - 1;
6747c6ae99SBarry Smith       for (ymin=0; ymin<dd->N; ymin++) {
6847c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
6947c6ae99SBarry Smith       }
7047c6ae99SBarry Smith     }
715b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
725b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
7347c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
7447c6ae99SBarry Smith 
755b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
7647c6ae99SBarry Smith     /* draw my box */
775b399a63SLisandro Dalcin     xmin = dd->xs/dd->w; xmax =(dd->xe-1)/dd->w; ymin = dd->ys; ymax = dd->ye - 1;
7847c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
7947c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
8047c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
8147c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
8247c6ae99SBarry Smith     /* put in numbers */
8347c6ae99SBarry Smith     base = (dd->base)/dd->w;
8447c6ae99SBarry Smith     for (y=ymin; y<=ymax; y++) {
8547c6ae99SBarry Smith       for (x=xmin; x<=xmax; x++) {
865b399a63SLisandro Dalcin         ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)base++);CHKERRQ(ierr);
8747c6ae99SBarry Smith         ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
8847c6ae99SBarry Smith       }
8947c6ae99SBarry Smith     }
905b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
915b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
9247c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
9347c6ae99SBarry Smith 
945b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
955b399a63SLisandro Dalcin     /* overlay ghost numbers, useful for error checking */
9645b6f7e9SBarry Smith     ierr = ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr);
975b399a63SLisandro Dalcin     base = 0; xmin = dd->Xs; xmax = dd->Xe; ymin = dd->Ys; ymax = dd->Ye;
9847c6ae99SBarry Smith     for (y=ymin; y<ymax; y++) {
9947c6ae99SBarry Smith       for (x=xmin; x<xmax; x++) {
10047c6ae99SBarry Smith         if ((base % dd->w) == 0) {
1015b399a63SLisandro Dalcin           ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)(idx[base/dd->w]));CHKERRQ(ierr);
10247c6ae99SBarry Smith           ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
10347c6ae99SBarry Smith         }
10447c6ae99SBarry Smith         base++;
10547c6ae99SBarry Smith       }
10647c6ae99SBarry Smith     }
107302440fdSBarry Smith     ierr = ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr);
1085b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
1095b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
11047c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
111*832b7cebSLisandro Dalcin     ierr = PetscDrawSave(draw);CHKERRQ(ierr);
1129a42bb27SBarry Smith   } else if (isbinary) {
1139a42bb27SBarry Smith     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
1149a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1159a42bb27SBarry Smith   } else if (ismatlab) {
1169a42bb27SBarry Smith     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
1179a42bb27SBarry Smith #endif
11811aeaf0aSBarry Smith   }
11947c6ae99SBarry Smith   PetscFunctionReturn(0);
12047c6ae99SBarry Smith }
12147c6ae99SBarry Smith 
12247c6ae99SBarry Smith /*
12347c6ae99SBarry Smith       M is number of grid points
12447c6ae99SBarry Smith       m is number of processors
12547c6ae99SBarry Smith 
12647c6ae99SBarry Smith */
12747c6ae99SBarry Smith #undef __FUNCT__
128aa219208SBarry Smith #define __FUNCT__ "DMDASplitComm2d"
1297087cfbeSBarry Smith PetscErrorCode  DMDASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm)
13047c6ae99SBarry Smith {
13147c6ae99SBarry Smith   PetscErrorCode ierr;
13247c6ae99SBarry Smith   PetscInt       m,n = 0,x = 0,y = 0;
13347c6ae99SBarry Smith   PetscMPIInt    size,csize,rank;
13447c6ae99SBarry Smith 
13547c6ae99SBarry Smith   PetscFunctionBegin;
13647c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
13747c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13847c6ae99SBarry Smith 
13947c6ae99SBarry Smith   csize = 4*size;
14047c6ae99SBarry Smith   do {
14147c6ae99SBarry 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);
14247c6ae99SBarry Smith     csize = csize/4;
14347c6ae99SBarry Smith 
144369cc0aeSBarry Smith     m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)csize)/((PetscReal)N)));
14547c6ae99SBarry Smith     if (!m) m = 1;
14647c6ae99SBarry Smith     while (m > 0) {
14747c6ae99SBarry Smith       n = csize/m;
14847c6ae99SBarry Smith       if (m*n == csize) break;
14947c6ae99SBarry Smith       m--;
15047c6ae99SBarry Smith     }
15147c6ae99SBarry Smith     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
15247c6ae99SBarry Smith 
15347c6ae99SBarry Smith     x = M/m + ((M % m) > ((csize-1) % m));
15447c6ae99SBarry Smith     y = (N + (csize-1)/m)/n;
15547c6ae99SBarry Smith   } while ((x < 4 || y < 4) && csize > 1);
15647c6ae99SBarry Smith   if (size != csize) {
15747c6ae99SBarry Smith     MPI_Group   entire_group,sub_group;
15847c6ae99SBarry Smith     PetscMPIInt i,*groupies;
15947c6ae99SBarry Smith 
16047c6ae99SBarry Smith     ierr = MPI_Comm_group(comm,&entire_group);CHKERRQ(ierr);
161785e854fSJed Brown     ierr = PetscMalloc1(csize,&groupies);CHKERRQ(ierr);
16247c6ae99SBarry Smith     for (i=0; i<csize; i++) {
16347c6ae99SBarry Smith       groupies[i] = (rank/csize)*csize + i;
16447c6ae99SBarry Smith     }
16547c6ae99SBarry Smith     ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);CHKERRQ(ierr);
16647c6ae99SBarry Smith     ierr = PetscFree(groupies);CHKERRQ(ierr);
16747c6ae99SBarry Smith     ierr = MPI_Comm_create(comm,sub_group,outcomm);CHKERRQ(ierr);
16847c6ae99SBarry Smith     ierr = MPI_Group_free(&entire_group);CHKERRQ(ierr);
16947c6ae99SBarry Smith     ierr = MPI_Group_free(&sub_group);CHKERRQ(ierr);
170aa219208SBarry Smith     ierr = PetscInfo1(0,"DMDASplitComm2d:Creating redundant coarse problems of size %d\n",csize);CHKERRQ(ierr);
17147c6ae99SBarry Smith   } else {
17247c6ae99SBarry Smith     *outcomm = comm;
17347c6ae99SBarry Smith   }
17447c6ae99SBarry Smith   PetscFunctionReturn(0);
17547c6ae99SBarry Smith }
17647c6ae99SBarry Smith 
17747c6ae99SBarry Smith #if defined(new)
17847c6ae99SBarry Smith #undef __FUNCT__
179aa219208SBarry Smith #define __FUNCT__ "DMDAGetDiagonal_MFFD"
18047c6ae99SBarry Smith /*
181aa219208SBarry Smith   DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
182aa219208SBarry Smith     function lives on a DMDA
18347c6ae99SBarry Smith 
18447c6ae99SBarry Smith         y ~= (F(u + ha) - F(u))/h,
18547c6ae99SBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
18647c6ae99SBarry Smith         u = current iterate
18747c6ae99SBarry Smith         h = difference interval
18847c6ae99SBarry Smith */
189aa219208SBarry Smith PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a)
19047c6ae99SBarry Smith {
19147c6ae99SBarry Smith   PetscScalar    h,*aa,*ww,v;
19247c6ae99SBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
19347c6ae99SBarry Smith   PetscErrorCode ierr;
19447c6ae99SBarry Smith   PetscInt       gI,nI;
19547c6ae99SBarry Smith   MatStencil     stencil;
196aa219208SBarry Smith   DMDALocalInfo  info;
19747c6ae99SBarry Smith 
19847c6ae99SBarry Smith   PetscFunctionBegin;
19947c6ae99SBarry Smith   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
20047c6ae99SBarry Smith   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
20147c6ae99SBarry Smith 
20247c6ae99SBarry Smith   ierr = VecGetArray(U,&ww);CHKERRQ(ierr);
20347c6ae99SBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
20447c6ae99SBarry Smith 
20547c6ae99SBarry Smith   nI = 0;
20647c6ae99SBarry Smith   h  = ww[gI];
20747c6ae99SBarry Smith   if (h == 0.0) h = 1.0;
20847c6ae99SBarry Smith   if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin;
20947c6ae99SBarry Smith   else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
21047c6ae99SBarry Smith   h *= epsilon;
21147c6ae99SBarry Smith 
21247c6ae99SBarry Smith   ww[gI] += h;
21347c6ae99SBarry Smith   ierr    = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
21447c6ae99SBarry Smith   aa[nI]  = (v - aa[nI])/h;
21547c6ae99SBarry Smith   ww[gI] -= h;
21647c6ae99SBarry Smith   nI++;
2178865f1eaSKarl Rupp 
21847c6ae99SBarry Smith   ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr);
21947c6ae99SBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
22047c6ae99SBarry Smith   PetscFunctionReturn(0);
22147c6ae99SBarry Smith }
22247c6ae99SBarry Smith #endif
22347c6ae99SBarry Smith 
22447c6ae99SBarry Smith #undef __FUNCT__
2259a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_2D"
2267087cfbeSBarry Smith PetscErrorCode  DMSetUp_DA_2D(DM da)
22747c6ae99SBarry Smith {
22847c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
22947c6ae99SBarry Smith   const PetscInt   M            = dd->M;
23047c6ae99SBarry Smith   const PetscInt   N            = dd->N;
23147c6ae99SBarry Smith   PetscInt         m            = dd->m;
23247c6ae99SBarry Smith   PetscInt         n            = dd->n;
23347c6ae99SBarry Smith   const PetscInt   dof          = dd->w;
23447c6ae99SBarry Smith   const PetscInt   s            = dd->s;
235bff4a2f0SMatthew G. Knepley   DMBoundaryType   bx           = dd->bx;
236bff4a2f0SMatthew G. Knepley   DMBoundaryType   by           = dd->by;
23719fd82e9SBarry Smith   DMDAStencilType  stencil_type = dd->stencil_type;
23847c6ae99SBarry Smith   PetscInt         *lx          = dd->lx;
23947c6ae99SBarry Smith   PetscInt         *ly          = dd->ly;
24047c6ae99SBarry Smith   MPI_Comm         comm;
24147c6ae99SBarry Smith   PetscMPIInt      rank,size;
242bd1fc5aeSBarry Smith   PetscInt         xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,IXs,IXe,IYs,IYe;
2438ea3bf28SBarry Smith   PetscInt         up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn;
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;
248bd1fc5aeSBarry Smith   VecScatter       gtol;
24945b6f7e9SBarry Smith   IS               to,from;
25047c6ae99SBarry Smith   PetscErrorCode   ierr;
25147c6ae99SBarry Smith 
25247c6ae99SBarry Smith   PetscFunctionBegin;
253bff4a2f0SMatthew G. Knepley   if (stencil_type == DMDA_STENCIL_BOX && (bx == DM_BOUNDARY_MIRROR || by == DM_BOUNDARY_MIRROR)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Mirror boundary and box stencil");
25447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
2553855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES)
2563855c12bSBarry 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);
2573855c12bSBarry Smith #endif
25847c6ae99SBarry Smith 
25947c6ae99SBarry Smith   if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
26047c6ae99SBarry Smith   if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
26147c6ae99SBarry Smith 
26247c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
26347c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
26447c6ae99SBarry Smith 
2657d310018SBarry Smith   dd->p = 1;
26647c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
26747c6ae99SBarry Smith     if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
26847c6ae99SBarry Smith     else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
26947c6ae99SBarry Smith   }
27047c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
27147c6ae99SBarry Smith     if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
27247c6ae99SBarry Smith     else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
27347c6ae99SBarry Smith   }
27447c6ae99SBarry Smith 
27547c6ae99SBarry Smith   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
27647c6ae99SBarry Smith     if (n != PETSC_DECIDE) {
27747c6ae99SBarry Smith       m = size/n;
27847c6ae99SBarry Smith     } else if (m != PETSC_DECIDE) {
27947c6ae99SBarry Smith       n = size/m;
28047c6ae99SBarry Smith     } else {
28147c6ae99SBarry Smith       /* try for squarish distribution */
282369cc0aeSBarry Smith       m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N)));
28347c6ae99SBarry Smith       if (!m) m = 1;
28447c6ae99SBarry Smith       while (m > 0) {
28547c6ae99SBarry Smith         n = size/m;
28647c6ae99SBarry Smith         if (m*n == size) break;
28747c6ae99SBarry Smith         m--;
28847c6ae99SBarry Smith       }
28947c6ae99SBarry Smith       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
29047c6ae99SBarry Smith     }
29147c6ae99SBarry 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 ");
29247c6ae99SBarry Smith   } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
29347c6ae99SBarry Smith 
29447c6ae99SBarry Smith   if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
29547c6ae99SBarry Smith   if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
29647c6ae99SBarry Smith 
29747c6ae99SBarry Smith   /*
29847c6ae99SBarry Smith      Determine locally owned region
29947c6ae99SBarry Smith      xs is the first local node number, x is the number of local nodes
30047c6ae99SBarry Smith   */
30147c6ae99SBarry Smith   if (!lx) {
302785e854fSJed Brown     ierr = PetscMalloc1(m, &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);
30647c6ae99SBarry Smith     }
30747c6ae99SBarry Smith   }
30847c6ae99SBarry Smith   x  = lx[rank % m];
30947c6ae99SBarry Smith   xs = 0;
31047c6ae99SBarry Smith   for (i=0; i<(rank % m); i++) {
31147c6ae99SBarry Smith     xs += lx[i];
31247c6ae99SBarry Smith   }
31347c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
31447c6ae99SBarry Smith   left = xs;
31547c6ae99SBarry Smith   for (i=(rank % m); i<m; i++) {
31647c6ae99SBarry Smith     left += lx[i];
31747c6ae99SBarry Smith   }
31847c6ae99SBarry 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);
31947c6ae99SBarry Smith #endif
32047c6ae99SBarry Smith 
32147c6ae99SBarry Smith   /*
32247c6ae99SBarry Smith      Determine locally owned region
32347c6ae99SBarry Smith      ys is the first local node number, y is the number of local nodes
32447c6ae99SBarry Smith   */
32547c6ae99SBarry Smith   if (!ly) {
326785e854fSJed Brown     ierr = PetscMalloc1(n, &dd->ly);CHKERRQ(ierr);
32747c6ae99SBarry Smith     ly   = dd->ly;
32847c6ae99SBarry Smith     for (i=0; i<n; i++) {
32947c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > i);
33047c6ae99SBarry Smith     }
33147c6ae99SBarry Smith   }
33247c6ae99SBarry Smith   y  = ly[rank/m];
33347c6ae99SBarry Smith   ys = 0;
33447c6ae99SBarry Smith   for (i=0; i<(rank/m); i++) {
33547c6ae99SBarry Smith     ys += ly[i];
33647c6ae99SBarry Smith   }
33747c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
33847c6ae99SBarry Smith   left = ys;
33947c6ae99SBarry Smith   for (i=(rank/m); i<n; i++) {
34047c6ae99SBarry Smith     left += ly[i];
34147c6ae99SBarry Smith   }
34247c6ae99SBarry 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);
34347c6ae99SBarry Smith #endif
34447c6ae99SBarry Smith 
345bcea557cSEthan Coon   /*
346bcea557cSEthan Coon    check if the scatter requires more than one process neighbor or wraps around
347bcea557cSEthan Coon    the domain more than once
348bcea557cSEthan Coon   */
349bff4a2f0SMatthew G. Knepley   if ((x < s) && ((m > 1) || (bx == DM_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s);
350bff4a2f0SMatthew G. Knepley   if ((y < s) && ((n > 1) || (by == DM_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s);
35147c6ae99SBarry Smith   xe = xs + x;
35247c6ae99SBarry Smith   ye = ys + y;
35347c6ae99SBarry Smith 
354ce00eea3SSatish Balay   /* determine ghost region (Xs) and region scattered into (IXs)  */
355d9c9ebe5SPeter Brune   if (xs-s > 0) {
356d9c9ebe5SPeter Brune     Xs = xs - s; IXs = xs - s;
35788661749SPeter Brune   } else {
35888661749SPeter Brune     if (bx) {
35988661749SPeter Brune       Xs = xs - s;
36088661749SPeter Brune     } else {
36188661749SPeter Brune       Xs = 0;
36288661749SPeter Brune     }
36388661749SPeter Brune     IXs = 0;
36488661749SPeter Brune   }
365d9c9ebe5SPeter Brune   if (xe+s <= M) {
366d9c9ebe5SPeter Brune     Xe = xe + s; IXe = xe + s;
36788661749SPeter Brune   } else {
36888661749SPeter Brune     if (bx) {
369d9c9ebe5SPeter Brune       Xs = xs - s; Xe = xe + s;
37088661749SPeter Brune     } else {
37188661749SPeter Brune       Xe = M;
37288661749SPeter Brune     }
37388661749SPeter Brune     IXe = M;
37488661749SPeter Brune   }
37547c6ae99SBarry Smith 
376bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) {
377d9c9ebe5SPeter Brune     IXs = xs - s;
378d9c9ebe5SPeter Brune     IXe = xe + s;
379d9c9ebe5SPeter Brune     Xs  = xs - s;
380d9c9ebe5SPeter Brune     Xe  = xe + s;
38188661749SPeter Brune   }
38247c6ae99SBarry Smith 
383d9c9ebe5SPeter Brune   if (ys-s > 0) {
384d9c9ebe5SPeter Brune     Ys = ys - s; IYs = ys - s;
38588661749SPeter Brune   } else {
38688661749SPeter Brune     if (by) {
38788661749SPeter Brune       Ys = ys - s;
38888661749SPeter Brune     } else {
38988661749SPeter Brune       Ys = 0;
39088661749SPeter Brune     }
39188661749SPeter Brune     IYs = 0;
39288661749SPeter Brune   }
393d9c9ebe5SPeter Brune   if (ye+s <= N) {
394d9c9ebe5SPeter Brune     Ye = ye + s; IYe = ye + s;
39588661749SPeter Brune   } else {
39688661749SPeter Brune     if (by) {
39788661749SPeter Brune       Ye = ye + s;
39888661749SPeter Brune     } else {
39988661749SPeter Brune       Ye = N;
40088661749SPeter Brune     }
40188661749SPeter Brune     IYe = N;
40288661749SPeter Brune   }
40388661749SPeter Brune 
404bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) {
405d9c9ebe5SPeter Brune     IYs = ys - s;
406d9c9ebe5SPeter Brune     IYe = ye + s;
407d9c9ebe5SPeter Brune     Ys  = ys - s;
408d9c9ebe5SPeter Brune     Ye  = ye + s;
40988661749SPeter Brune   }
41088661749SPeter Brune 
41188661749SPeter Brune   /* stencil length in each direction */
412d9c9ebe5SPeter Brune   s_x = s;
413d9c9ebe5SPeter Brune   s_y = s;
41447c6ae99SBarry Smith 
41547c6ae99SBarry Smith   /* determine starting point of each processor */
41647c6ae99SBarry Smith   nn       = x*y;
417dcca6d9dSJed Brown   ierr     = PetscMalloc2(size+1,&bases,size,&ldims);CHKERRQ(ierr);
41847c6ae99SBarry Smith   ierr     = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
41947c6ae99SBarry Smith   bases[0] = 0;
42047c6ae99SBarry Smith   for (i=1; i<=size; i++) {
42147c6ae99SBarry Smith     bases[i] = ldims[i-1];
42247c6ae99SBarry Smith   }
42347c6ae99SBarry Smith   for (i=1; i<=size; i++) {
42447c6ae99SBarry Smith     bases[i] += bases[i-1];
42547c6ae99SBarry Smith   }
426ce00eea3SSatish Balay   base = bases[rank]*dof;
42747c6ae99SBarry Smith 
42847c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
429ce00eea3SSatish Balay   dd->Nlocal = x*y*dof;
430b1fb7eb7SBarry Smith   ierr       = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global);CHKERRQ(ierr);
431ce00eea3SSatish Balay   dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof;
432b1fb7eb7SBarry Smith   ierr       = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local);CHKERRQ(ierr);
43347c6ae99SBarry Smith 
43447c6ae99SBarry Smith   /* generate appropriate vector scatters */
43547c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
43602fe608eSBarry Smith   ierr  = PetscMalloc1((IXe-IXs)*(IYe-IYs),&idx);CHKERRQ(ierr);
437ce00eea3SSatish Balay   left  = xs - Xs; right = left + x;
438ce00eea3SSatish Balay   down  = ys - Ys; up = down + y;
43947c6ae99SBarry Smith   count = 0;
44047c6ae99SBarry Smith   for (i=down; i<up; i++) {
441ce00eea3SSatish Balay     for (j=left; j<right; j++) {
442ce00eea3SSatish Balay       idx[count++] = i*(Xe-Xs) + j;
44347c6ae99SBarry Smith     }
44447c6ae99SBarry Smith   }
44547c6ae99SBarry Smith 
446ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
447ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
448d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_BOX) {
449ce00eea3SSatish Balay     left  = IXs - Xs; right = left + (IXe-IXs);
450ce00eea3SSatish Balay     down  = IYs - Ys; up = down + (IYe-IYs);
451ce00eea3SSatish Balay     count = 0;
452ce00eea3SSatish Balay     for (i=down; i<up; i++) {
453ce00eea3SSatish Balay       for (j=left; j<right; j++) {
454ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
455ce00eea3SSatish Balay       }
456ce00eea3SSatish Balay     }
457ce00eea3SSatish Balay     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
458ce00eea3SSatish Balay 
45947c6ae99SBarry Smith   } else {
46047c6ae99SBarry Smith     /* must drop into cross shape region */
46147c6ae99SBarry Smith     /*       ---------|
46247c6ae99SBarry Smith             |  top    |
463ce00eea3SSatish Balay          |---         ---| up
46447c6ae99SBarry Smith          |   middle      |
46547c6ae99SBarry Smith          |               |
466ce00eea3SSatish Balay          ----         ---- down
46747c6ae99SBarry Smith             | bottom  |
46847c6ae99SBarry Smith             -----------
46947c6ae99SBarry Smith          Xs xs        xe Xe */
470ce00eea3SSatish Balay     left  = xs - Xs; right = left + x;
471ce00eea3SSatish Balay     down  = ys - Ys; up = down + y;
47247c6ae99SBarry Smith     count = 0;
473ce00eea3SSatish Balay     /* bottom */
474ce00eea3SSatish Balay     for (i=(IYs-Ys); i<down; i++) {
475ce00eea3SSatish Balay       for (j=left; j<right; j++) {
476ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
47747c6ae99SBarry Smith       }
47847c6ae99SBarry Smith     }
47947c6ae99SBarry Smith     /* middle */
48047c6ae99SBarry Smith     for (i=down; i<up; i++) {
481ce00eea3SSatish Balay       for (j=(IXs-Xs); j<(IXe-Xs); j++) {
482ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
48347c6ae99SBarry Smith       }
48447c6ae99SBarry Smith     }
48547c6ae99SBarry Smith     /* top */
486ce00eea3SSatish Balay     for (i=up; i<up+IYe-ye; i++) {
487ce00eea3SSatish Balay       for (j=left; j<right; j++) {
488ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
48947c6ae99SBarry Smith       }
49047c6ae99SBarry Smith     }
49147c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
49247c6ae99SBarry Smith   }
49347c6ae99SBarry Smith 
49447c6ae99SBarry Smith 
49547c6ae99SBarry Smith   /* determine who lies on each side of us stored in    n6 n7 n8
49647c6ae99SBarry Smith                                                         n3    n5
49747c6ae99SBarry Smith                                                         n0 n1 n2
49847c6ae99SBarry Smith   */
49947c6ae99SBarry Smith 
50047c6ae99SBarry Smith   /* Assume the Non-Periodic Case */
50147c6ae99SBarry Smith   n1 = rank - m;
50247c6ae99SBarry Smith   if (rank % m) {
50347c6ae99SBarry Smith     n0 = n1 - 1;
50447c6ae99SBarry Smith   } else {
50547c6ae99SBarry Smith     n0 = -1;
50647c6ae99SBarry Smith   }
50747c6ae99SBarry Smith   if ((rank+1) % m) {
50847c6ae99SBarry Smith     n2 = n1 + 1;
50947c6ae99SBarry Smith     n5 = rank + 1;
51047c6ae99SBarry Smith     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
51147c6ae99SBarry Smith   } else {
51247c6ae99SBarry Smith     n2 = -1; n5 = -1; n8 = -1;
51347c6ae99SBarry Smith   }
51447c6ae99SBarry Smith   if (rank % m) {
51547c6ae99SBarry Smith     n3 = rank - 1;
51647c6ae99SBarry Smith     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
51747c6ae99SBarry Smith   } else {
51847c6ae99SBarry Smith     n3 = -1; n6 = -1;
51947c6ae99SBarry Smith   }
52047c6ae99SBarry Smith   n7 = rank + m; if (n7 >= m*n) n7 = -1;
52147c6ae99SBarry Smith 
522bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC && by == DM_BOUNDARY_PERIODIC) {
52347c6ae99SBarry Smith     /* Modify for Periodic Cases */
52447c6ae99SBarry Smith     /* Handle all four corners */
52547c6ae99SBarry Smith     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
52647c6ae99SBarry Smith     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
52747c6ae99SBarry Smith     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
52847c6ae99SBarry Smith     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
52947c6ae99SBarry Smith 
53047c6ae99SBarry Smith     /* Handle Top and Bottom Sides */
53147c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
53247c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
53347c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
53447c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
53547c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
53647c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
53747c6ae99SBarry Smith 
53847c6ae99SBarry Smith     /* Handle Left and Right Sides */
53947c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
54047c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
54147c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
54247c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
54347c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
54447c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
545bff4a2f0SMatthew G. Knepley   } else if (by == DM_BOUNDARY_PERIODIC) {  /* Handle Top and Bottom Sides */
546ce00eea3SSatish Balay     if (n1 < 0) n1 = rank + m * (n-1);
547ce00eea3SSatish Balay     if (n7 < 0) n7 = rank - m * (n-1);
548ce00eea3SSatish Balay     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
549ce00eea3SSatish Balay     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
550ce00eea3SSatish Balay     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
551ce00eea3SSatish Balay     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
552bff4a2f0SMatthew G. Knepley   } else if (bx == DM_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */
553ce00eea3SSatish Balay     if (n3 < 0) n3 = rank + (m-1);
554ce00eea3SSatish Balay     if (n5 < 0) n5 = rank - (m-1);
555ce00eea3SSatish Balay     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
556ce00eea3SSatish Balay     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
557ce00eea3SSatish Balay     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
558ce00eea3SSatish Balay     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
55947c6ae99SBarry Smith   }
560ce00eea3SSatish Balay 
561785e854fSJed Brown   ierr = PetscMalloc1(9,&dd->neighbors);CHKERRQ(ierr);
5628865f1eaSKarl Rupp 
56347c6ae99SBarry Smith   dd->neighbors[0] = n0;
56447c6ae99SBarry Smith   dd->neighbors[1] = n1;
56547c6ae99SBarry Smith   dd->neighbors[2] = n2;
56647c6ae99SBarry Smith   dd->neighbors[3] = n3;
56747c6ae99SBarry Smith   dd->neighbors[4] = rank;
56847c6ae99SBarry Smith   dd->neighbors[5] = n5;
56947c6ae99SBarry Smith   dd->neighbors[6] = n6;
57047c6ae99SBarry Smith   dd->neighbors[7] = n7;
57147c6ae99SBarry Smith   dd->neighbors[8] = n8;
57247c6ae99SBarry Smith 
573d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
57447c6ae99SBarry Smith     /* save corner processor numbers */
57547c6ae99SBarry Smith     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
57647c6ae99SBarry Smith     n0  = n2 = n6 = n8 = -1;
57747c6ae99SBarry Smith   }
57847c6ae99SBarry Smith 
579785e854fSJed Brown   ierr = PetscMalloc1((Xe-Xs)*(Ye-Ys),&idx);CHKERRQ(ierr);
58047c6ae99SBarry Smith 
581ce00eea3SSatish Balay   nn = 0;
58247c6ae99SBarry Smith   xbase = bases[rank];
58347c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
58447c6ae99SBarry Smith     if (n0 >= 0) { /* left below */
585ce00eea3SSatish Balay       x_t = lx[n0 % m];
58647c6ae99SBarry Smith       y_t = ly[(n0/m)];
58747c6ae99SBarry Smith       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
5888865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
58947c6ae99SBarry Smith     }
590ac119b13SBarry Smith 
59147c6ae99SBarry Smith     if (n1 >= 0) { /* directly below */
59247c6ae99SBarry Smith       x_t = x;
59347c6ae99SBarry Smith       y_t = ly[(n1/m)];
59447c6ae99SBarry Smith       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
5958865f1eaSKarl Rupp       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
596bff4a2f0SMatthew G. Knepley     } else if (by == DM_BOUNDARY_MIRROR) {
5978865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
59847c6ae99SBarry Smith     }
599ac119b13SBarry Smith 
60047c6ae99SBarry Smith     if (n2 >= 0) { /* right below */
601ce00eea3SSatish Balay       x_t = lx[n2 % m];
60247c6ae99SBarry Smith       y_t = ly[(n2/m)];
60347c6ae99SBarry Smith       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
6048865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
60547c6ae99SBarry Smith     }
60647c6ae99SBarry Smith   }
60747c6ae99SBarry Smith 
60847c6ae99SBarry Smith   for (i=0; i<y; i++) {
60947c6ae99SBarry Smith     if (n3 >= 0) { /* directly left */
610ce00eea3SSatish Balay       x_t = lx[n3 % m];
61147c6ae99SBarry Smith       /* y_t = y; */
61247c6ae99SBarry Smith       s_t = bases[n3] + (i+1)*x_t - s_x;
6138865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
614bff4a2f0SMatthew G. Knepley     } else if (bx == DM_BOUNDARY_MIRROR) {
6158865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
61647c6ae99SBarry Smith     }
61747c6ae99SBarry Smith 
6188865f1eaSKarl Rupp     for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
61947c6ae99SBarry Smith 
62047c6ae99SBarry Smith     if (n5 >= 0) { /* directly right */
621ce00eea3SSatish Balay       x_t = lx[n5 % m];
62247c6ae99SBarry Smith       /* y_t = y; */
62347c6ae99SBarry Smith       s_t = bases[n5] + (i)*x_t;
6248865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
625bff4a2f0SMatthew G. Knepley     } else if (bx == DM_BOUNDARY_MIRROR) {
6268865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
62747c6ae99SBarry Smith     }
62847c6ae99SBarry Smith   }
62947c6ae99SBarry Smith 
63047c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
63147c6ae99SBarry Smith     if (n6 >= 0) { /* left above */
632ce00eea3SSatish Balay       x_t = lx[n6 % m];
63347c6ae99SBarry Smith       /* y_t = ly[(n6/m)]; */
63447c6ae99SBarry Smith       s_t = bases[n6] + (i)*x_t - s_x;
6358865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
63647c6ae99SBarry Smith     }
637ac119b13SBarry Smith 
63847c6ae99SBarry Smith     if (n7 >= 0) { /* directly above */
63947c6ae99SBarry Smith       x_t = x;
64047c6ae99SBarry Smith       /* y_t = ly[(n7/m)]; */
64147c6ae99SBarry Smith       s_t = bases[n7] + (i-1)*x_t;
6428865f1eaSKarl Rupp       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
643bff4a2f0SMatthew G. Knepley     } else if (by == DM_BOUNDARY_MIRROR) {
6448865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
64547c6ae99SBarry Smith     }
646ac119b13SBarry Smith 
64747c6ae99SBarry Smith     if (n8 >= 0) { /* right above */
648ce00eea3SSatish Balay       x_t = lx[n8 % m];
64947c6ae99SBarry Smith       /* y_t = ly[(n8/m)]; */
65047c6ae99SBarry Smith       s_t = bases[n8] + (i-1)*x_t;
6518865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
65247c6ae99SBarry Smith     }
65347c6ae99SBarry Smith   }
65447c6ae99SBarry Smith 
655b1fb7eb7SBarry Smith   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from);CHKERRQ(ierr);
65647c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
6573bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr);
658fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
659fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
66047c6ae99SBarry Smith 
661d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
662ce00eea3SSatish Balay     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
663ce00eea3SSatish Balay   }
664ce00eea3SSatish Balay 
66588661749SPeter Brune   if (((stencil_type == DMDA_STENCIL_STAR)  ||
666bff4a2f0SMatthew G. Knepley        (bx && bx != DM_BOUNDARY_PERIODIC) ||
667bff4a2f0SMatthew G. Knepley        (by && by != DM_BOUNDARY_PERIODIC))) {
66847c6ae99SBarry Smith     /*
66947c6ae99SBarry Smith         Recompute the local to global mappings, this time keeping the
670ce00eea3SSatish Balay       information about the cross corner processor numbers and any ghosted
671ce00eea3SSatish Balay       but not periodic indices.
67247c6ae99SBarry Smith     */
67347c6ae99SBarry Smith     nn    = 0;
67447c6ae99SBarry Smith     xbase = bases[rank];
67547c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
67647c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
677ce00eea3SSatish Balay         x_t = lx[n0 % m];
67847c6ae99SBarry Smith         y_t = ly[(n0/m)];
67947c6ae99SBarry Smith         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
6808865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
681ce00eea3SSatish Balay       } else if (xs-Xs > 0 && ys-Ys > 0) {
6828865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
68347c6ae99SBarry Smith       }
68447c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
68547c6ae99SBarry Smith         x_t = x;
68647c6ae99SBarry Smith         y_t = ly[(n1/m)];
68747c6ae99SBarry Smith         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
6888865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
689ce00eea3SSatish Balay       } else if (ys-Ys > 0) {
690bff4a2f0SMatthew G. Knepley         if (by == DM_BOUNDARY_MIRROR) {
6918865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
692624904c4SBarry Smith         } else {
6938865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
69447c6ae99SBarry Smith         }
695624904c4SBarry Smith       }
69647c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
697ce00eea3SSatish Balay         x_t = lx[n2 % m];
69847c6ae99SBarry Smith         y_t = ly[(n2/m)];
69947c6ae99SBarry Smith         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
7008865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
701ce00eea3SSatish Balay       } else if (Xe-xe> 0 && ys-Ys > 0) {
7028865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
70347c6ae99SBarry Smith       }
70447c6ae99SBarry Smith     }
70547c6ae99SBarry Smith 
70647c6ae99SBarry Smith     for (i=0; i<y; i++) {
70747c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
708ce00eea3SSatish Balay         x_t = lx[n3 % m];
70947c6ae99SBarry Smith         /* y_t = y; */
71047c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x;
7118865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
712ce00eea3SSatish Balay       } else if (xs-Xs > 0) {
713bff4a2f0SMatthew G. Knepley         if (bx == DM_BOUNDARY_MIRROR) {
7148865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
715624904c4SBarry Smith         } else {
7168865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
71747c6ae99SBarry Smith         }
718624904c4SBarry Smith       }
71947c6ae99SBarry Smith 
7208865f1eaSKarl Rupp       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;
7268865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
727ce00eea3SSatish Balay       } else if (Xe-xe > 0) {
728bff4a2f0SMatthew G. Knepley         if (bx == DM_BOUNDARY_MIRROR) {
7298865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
730624904c4SBarry Smith         } else {
7318865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
73247c6ae99SBarry Smith         }
73347c6ae99SBarry Smith       }
734624904c4SBarry Smith     }
73547c6ae99SBarry Smith 
73647c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
73747c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
738ce00eea3SSatish Balay         x_t = lx[n6 % m];
73947c6ae99SBarry Smith         /* y_t = ly[(n6/m)]; */
74047c6ae99SBarry Smith         s_t = bases[n6] + (i)*x_t - s_x;
7418865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
742ce00eea3SSatish Balay       } else if (xs-Xs > 0 && Ye-ye > 0) {
7438865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
74447c6ae99SBarry Smith       }
74547c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
74647c6ae99SBarry Smith         x_t = x;
74747c6ae99SBarry Smith         /* y_t = ly[(n7/m)]; */
74847c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t;
7498865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
750ce00eea3SSatish Balay       } else if (Ye-ye > 0) {
751bff4a2f0SMatthew G. Knepley         if (by == DM_BOUNDARY_MIRROR) {
7528865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
753624904c4SBarry Smith         } else {
7548865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
75547c6ae99SBarry Smith         }
756624904c4SBarry Smith       }
75747c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
758ce00eea3SSatish Balay         x_t = lx[n8 % m];
75947c6ae99SBarry Smith         /* y_t = ly[(n8/m)]; */
76047c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t;
7618865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
762ce00eea3SSatish Balay       } else if (Xe-xe > 0 && Ye-ye > 0) {
7638865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
76447c6ae99SBarry Smith       }
76547c6ae99SBarry Smith     }
76647c6ae99SBarry Smith   }
767ce00eea3SSatish Balay   /*
768ce00eea3SSatish Balay      Set the local to global ordering in the global vector, this allows use
769ce00eea3SSatish Balay      of VecSetValuesLocal().
770ce00eea3SSatish Balay   */
77145b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr);
7723bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr);
77347c6ae99SBarry Smith 
774ce00eea3SSatish Balay   ierr  = PetscFree2(bases,ldims);CHKERRQ(ierr);
77547c6ae99SBarry Smith   dd->m = m;  dd->n  = n;
776ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
777ce00eea3SSatish Balay   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
778ce00eea3SSatish Balay   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
77947c6ae99SBarry Smith 
780fcfd50ebSBarry Smith   ierr = VecDestroy(&local);CHKERRQ(ierr);
781fcfd50ebSBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
78247c6ae99SBarry Smith 
78347c6ae99SBarry Smith   dd->gtol      = gtol;
78447c6ae99SBarry Smith   dd->base      = base;
7859a42bb27SBarry Smith   da->ops->view = DMView_DA_2d;
7860298fd71SBarry Smith   dd->ltol      = NULL;
7870298fd71SBarry Smith   dd->ao        = NULL;
78847c6ae99SBarry Smith   PetscFunctionReturn(0);
78947c6ae99SBarry Smith }
79047c6ae99SBarry Smith 
79147c6ae99SBarry Smith #undef __FUNCT__
792aa219208SBarry Smith #define __FUNCT__ "DMDACreate2d"
79347c6ae99SBarry Smith /*@C
794aa219208SBarry Smith    DMDACreate2d -  Creates an object that will manage the communication of  two-dimensional
79547c6ae99SBarry Smith    regular array data that is distributed across some processors.
79647c6ae99SBarry Smith 
79747c6ae99SBarry Smith    Collective on MPI_Comm
79847c6ae99SBarry Smith 
79947c6ae99SBarry Smith    Input Parameters:
80047c6ae99SBarry Smith +  comm - MPI communicator
8011321219cSEthan Coon .  bx,by - type of ghost nodes the array have.
802bff4a2f0SMatthew G. Knepley          Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC.
803aa219208SBarry Smith .  stencil_type - stencil type.  Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
80447c6ae99SBarry 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
80547c6ae99SBarry Smith             from the command line with -da_grid_x <M> -da_grid_y <N>)
80647c6ae99SBarry Smith .  m,n - corresponding number of processors in each dimension
80747c6ae99SBarry Smith          (or PETSC_DECIDE to have calculated)
80847c6ae99SBarry Smith .  dof - number of degrees of freedom per node
80947c6ae99SBarry Smith .  s - stencil width
81047c6ae99SBarry Smith -  lx, ly - arrays containing the number of nodes in each cell along
8110298fd71SBarry Smith            the x and y coordinates, or NULL. If non-null, these
81247c6ae99SBarry Smith            must be of length as m and n, and the corresponding
81347c6ae99SBarry Smith            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
81447c6ae99SBarry Smith            must be M, and the sum of the ly[] entries must be N.
81547c6ae99SBarry Smith 
81647c6ae99SBarry Smith    Output Parameter:
81747c6ae99SBarry Smith .  da - the resulting distributed array object
81847c6ae99SBarry Smith 
81947c6ae99SBarry Smith    Options Database Key:
820706a11cbSBarry Smith +  -dm_view - Calls DMView() at the conclusion of DMDACreate2d()
82147c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
82247c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
82347c6ae99SBarry Smith .  -da_processors_x <nx> - number of processors in x direction
82447c6ae99SBarry Smith .  -da_processors_y <ny> - number of processors in y direction
825e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement ratio in x direction
826e0f5d30fSBarry Smith .  -da_refine_y <ry> - refinement ratio in y direction
827e0f5d30fSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating, if M or N < 0
828e0f5d30fSBarry Smith 
82947c6ae99SBarry Smith 
83047c6ae99SBarry Smith    Level: beginner
83147c6ae99SBarry Smith 
83247c6ae99SBarry Smith    Notes:
833aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
834aa219208SBarry Smith    standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
83547c6ae99SBarry Smith    the standard 9-pt stencil.
83647c6ae99SBarry Smith 
837aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
838564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
839564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
84047c6ae99SBarry Smith 
84147c6ae99SBarry Smith .keywords: distributed array, create, two-dimensional
84247c6ae99SBarry Smith 
843aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
84499f0b487SRichard Tran Mills           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(),
845d461ba97SBarry Smith           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
84647c6ae99SBarry Smith 
84747c6ae99SBarry Smith @*/
848fe16a2e9SBarry Smith 
849bff4a2f0SMatthew G. Knepley PetscErrorCode  DMDACreate2d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMDAStencilType stencil_type,
8509a42bb27SBarry Smith                           PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
85147c6ae99SBarry Smith {
85247c6ae99SBarry Smith   PetscErrorCode ierr;
85347c6ae99SBarry Smith 
85447c6ae99SBarry Smith   PetscFunctionBegin;
855aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
856c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*da, 2);CHKERRQ(ierr);
857aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, 1);CHKERRQ(ierr);
858aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr);
859bff4a2f0SMatthew G. Knepley   ierr = DMDASetBoundaryType(*da, bx, by, DM_BOUNDARY_NONE);CHKERRQ(ierr);
860aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
861aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
862aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
8630298fd71SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, NULL);CHKERRQ(ierr);
86447c6ae99SBarry Smith   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
8659a42bb27SBarry Smith   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
8669a42bb27SBarry Smith   ierr = DMSetUp(*da);CHKERRQ(ierr);
86747c6ae99SBarry Smith   PetscFunctionReturn(0);
86847c6ae99SBarry Smith }
869