xref: /petsc/src/dm/impls/da/da2.c (revision 45b6f7e9ecd564e8bb535880b270d4a8084ff211)
19a42bb27SBarry Smith 
24035e84dSBarry 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);
337b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);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);
377b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);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);
5247c6ae99SBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
536636e97aSMatthew G Knepley     if (!da->coordinates) {
5447c6ae99SBarry Smith       ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
5547c6ae99SBarry Smith     }
5647c6ae99SBarry Smith     ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
5747c6ae99SBarry Smith 
5847c6ae99SBarry Smith     /* first processor draw all node lines */
5947c6ae99SBarry Smith     if (!rank) {
6047c6ae99SBarry Smith       ymin = 0.0; ymax = dd->N - 1;
6147c6ae99SBarry Smith       for (xmin=0; xmin<dd->M; xmin++) {
6247c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
6347c6ae99SBarry Smith       }
6447c6ae99SBarry Smith       xmin = 0.0; xmax = dd->M - 1;
6547c6ae99SBarry Smith       for (ymin=0; ymin<dd->N; ymin++) {
6647c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
6747c6ae99SBarry Smith       }
6847c6ae99SBarry Smith     }
6947c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
7047c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
7147c6ae99SBarry Smith 
7247c6ae99SBarry Smith     /* draw my box */
7347c6ae99SBarry Smith     ymin = dd->ys; ymax = dd->ye - 1; xmin = dd->xs/dd->w;
7447c6ae99SBarry Smith     xmax =(dd->xe-1)/dd->w;
7547c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
7647c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
7747c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
7847c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
7947c6ae99SBarry Smith 
8047c6ae99SBarry Smith     /* put in numbers */
8147c6ae99SBarry Smith     base = (dd->base)/dd->w;
8247c6ae99SBarry Smith     for (y=ymin; y<=ymax; y++) {
8347c6ae99SBarry Smith       for (x=xmin; x<=xmax; x++) {
8447c6ae99SBarry Smith         sprintf(node,"%d",(int)base++);
8547c6ae99SBarry Smith         ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
8647c6ae99SBarry Smith       }
8747c6ae99SBarry Smith     }
8847c6ae99SBarry Smith 
8947c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
9047c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
9147c6ae99SBarry Smith     /* overlay ghost numbers, useful for error checking */
9247c6ae99SBarry Smith     /* put in numbers */
9347c6ae99SBarry Smith 
948ea3bf28SBarry Smith     base = 0;
95*45b6f7e9SBarry Smith     ierr = ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr);
9647c6ae99SBarry Smith     ymin = dd->Ys; ymax = dd->Ye; xmin = dd->Xs; xmax = dd->Xe;
9747c6ae99SBarry Smith     for (y=ymin; y<ymax; y++) {
9847c6ae99SBarry Smith       for (x=xmin; x<xmax; x++) {
9947c6ae99SBarry Smith         if ((base % dd->w) == 0) {
100bfe97906SBarry Smith           sprintf(node,"%d",(int)(idx[base/dd->w]));
10147c6ae99SBarry Smith           ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
10247c6ae99SBarry Smith         }
10347c6ae99SBarry Smith         base++;
10447c6ae99SBarry Smith       }
10547c6ae99SBarry Smith     }
106*45b6f7e9SBarry Smith     ierr = ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx);
10747c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
10847c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
1099a42bb27SBarry Smith   } else if (isbinary) {
1109a42bb27SBarry Smith     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
1119a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1129a42bb27SBarry Smith   } else if (ismatlab) {
1139a42bb27SBarry Smith     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
1149a42bb27SBarry Smith #endif
11511aeaf0aSBarry Smith   }
11647c6ae99SBarry Smith   PetscFunctionReturn(0);
11747c6ae99SBarry Smith }
11847c6ae99SBarry Smith 
11947c6ae99SBarry Smith /*
12047c6ae99SBarry Smith       M is number of grid points
12147c6ae99SBarry Smith       m is number of processors
12247c6ae99SBarry Smith 
12347c6ae99SBarry Smith */
12447c6ae99SBarry Smith #undef __FUNCT__
125aa219208SBarry Smith #define __FUNCT__ "DMDASplitComm2d"
1267087cfbeSBarry Smith PetscErrorCode  DMDASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm)
12747c6ae99SBarry Smith {
12847c6ae99SBarry Smith   PetscErrorCode ierr;
12947c6ae99SBarry Smith   PetscInt       m,n = 0,x = 0,y = 0;
13047c6ae99SBarry Smith   PetscMPIInt    size,csize,rank;
13147c6ae99SBarry Smith 
13247c6ae99SBarry Smith   PetscFunctionBegin;
13347c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
13447c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13547c6ae99SBarry Smith 
13647c6ae99SBarry Smith   csize = 4*size;
13747c6ae99SBarry Smith   do {
13847c6ae99SBarry 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);
13947c6ae99SBarry Smith     csize = csize/4;
14047c6ae99SBarry Smith 
141369cc0aeSBarry Smith     m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)csize)/((PetscReal)N)));
14247c6ae99SBarry Smith     if (!m) m = 1;
14347c6ae99SBarry Smith     while (m > 0) {
14447c6ae99SBarry Smith       n = csize/m;
14547c6ae99SBarry Smith       if (m*n == csize) break;
14647c6ae99SBarry Smith       m--;
14747c6ae99SBarry Smith     }
14847c6ae99SBarry Smith     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
14947c6ae99SBarry Smith 
15047c6ae99SBarry Smith     x = M/m + ((M % m) > ((csize-1) % m));
15147c6ae99SBarry Smith     y = (N + (csize-1)/m)/n;
15247c6ae99SBarry Smith   } while ((x < 4 || y < 4) && csize > 1);
15347c6ae99SBarry Smith   if (size != csize) {
15447c6ae99SBarry Smith     MPI_Group   entire_group,sub_group;
15547c6ae99SBarry Smith     PetscMPIInt i,*groupies;
15647c6ae99SBarry Smith 
15747c6ae99SBarry Smith     ierr = MPI_Comm_group(comm,&entire_group);CHKERRQ(ierr);
158785e854fSJed Brown     ierr = PetscMalloc1(csize,&groupies);CHKERRQ(ierr);
15947c6ae99SBarry Smith     for (i=0; i<csize; i++) {
16047c6ae99SBarry Smith       groupies[i] = (rank/csize)*csize + i;
16147c6ae99SBarry Smith     }
16247c6ae99SBarry Smith     ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);CHKERRQ(ierr);
16347c6ae99SBarry Smith     ierr = PetscFree(groupies);CHKERRQ(ierr);
16447c6ae99SBarry Smith     ierr = MPI_Comm_create(comm,sub_group,outcomm);CHKERRQ(ierr);
16547c6ae99SBarry Smith     ierr = MPI_Group_free(&entire_group);CHKERRQ(ierr);
16647c6ae99SBarry Smith     ierr = MPI_Group_free(&sub_group);CHKERRQ(ierr);
167aa219208SBarry Smith     ierr = PetscInfo1(0,"DMDASplitComm2d:Creating redundant coarse problems of size %d\n",csize);CHKERRQ(ierr);
16847c6ae99SBarry Smith   } else {
16947c6ae99SBarry Smith     *outcomm = comm;
17047c6ae99SBarry Smith   }
17147c6ae99SBarry Smith   PetscFunctionReturn(0);
17247c6ae99SBarry Smith }
17347c6ae99SBarry Smith 
17447c6ae99SBarry Smith #if defined(new)
17547c6ae99SBarry Smith #undef __FUNCT__
176aa219208SBarry Smith #define __FUNCT__ "DMDAGetDiagonal_MFFD"
17747c6ae99SBarry Smith /*
178aa219208SBarry Smith   DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
179aa219208SBarry Smith     function lives on a DMDA
18047c6ae99SBarry Smith 
18147c6ae99SBarry Smith         y ~= (F(u + ha) - F(u))/h,
18247c6ae99SBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
18347c6ae99SBarry Smith         u = current iterate
18447c6ae99SBarry Smith         h = difference interval
18547c6ae99SBarry Smith */
186aa219208SBarry Smith PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a)
18747c6ae99SBarry Smith {
18847c6ae99SBarry Smith   PetscScalar    h,*aa,*ww,v;
18947c6ae99SBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
19047c6ae99SBarry Smith   PetscErrorCode ierr;
19147c6ae99SBarry Smith   PetscInt       gI,nI;
19247c6ae99SBarry Smith   MatStencil     stencil;
193aa219208SBarry Smith   DMDALocalInfo  info;
19447c6ae99SBarry Smith 
19547c6ae99SBarry Smith   PetscFunctionBegin;
19647c6ae99SBarry Smith   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
19747c6ae99SBarry Smith   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
19847c6ae99SBarry Smith 
19947c6ae99SBarry Smith   ierr = VecGetArray(U,&ww);CHKERRQ(ierr);
20047c6ae99SBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
20147c6ae99SBarry Smith 
20247c6ae99SBarry Smith   nI = 0;
20347c6ae99SBarry Smith   h  = ww[gI];
20447c6ae99SBarry Smith   if (h == 0.0) h = 1.0;
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   h *= epsilon;
20847c6ae99SBarry Smith 
20947c6ae99SBarry Smith   ww[gI] += h;
21047c6ae99SBarry Smith   ierr    = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
21147c6ae99SBarry Smith   aa[nI]  = (v - aa[nI])/h;
21247c6ae99SBarry Smith   ww[gI] -= h;
21347c6ae99SBarry Smith   nI++;
2148865f1eaSKarl Rupp 
21547c6ae99SBarry Smith   ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr);
21647c6ae99SBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
21747c6ae99SBarry Smith   PetscFunctionReturn(0);
21847c6ae99SBarry Smith }
21947c6ae99SBarry Smith #endif
22047c6ae99SBarry Smith 
22147c6ae99SBarry Smith #undef __FUNCT__
2229a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_2D"
2237087cfbeSBarry Smith PetscErrorCode  DMSetUp_DA_2D(DM da)
22447c6ae99SBarry Smith {
22547c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
22647c6ae99SBarry Smith   const PetscInt   M            = dd->M;
22747c6ae99SBarry Smith   const PetscInt   N            = dd->N;
22847c6ae99SBarry Smith   PetscInt         m            = dd->m;
22947c6ae99SBarry Smith   PetscInt         n            = dd->n;
23047c6ae99SBarry Smith   const PetscInt   dof          = dd->w;
23147c6ae99SBarry Smith   const PetscInt   s            = dd->s;
232bff4a2f0SMatthew G. Knepley   DMBoundaryType   bx           = dd->bx;
233bff4a2f0SMatthew G. Knepley   DMBoundaryType   by           = dd->by;
23419fd82e9SBarry Smith   DMDAStencilType  stencil_type = dd->stencil_type;
23547c6ae99SBarry Smith   PetscInt         *lx          = dd->lx;
23647c6ae99SBarry Smith   PetscInt         *ly          = dd->ly;
23747c6ae99SBarry Smith   MPI_Comm         comm;
23847c6ae99SBarry Smith   PetscMPIInt      rank,size;
239ce00eea3SSatish Balay   PetscInt         xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,start,end,IXs,IXe,IYs,IYe;
2408ea3bf28SBarry Smith   PetscInt         up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn;
241db87c5ecSEthan Coon   PetscInt         xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
24247c6ae99SBarry Smith   PetscInt         s_x,s_y; /* s proportionalized to w */
24347c6ae99SBarry Smith   PetscInt         sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
24447c6ae99SBarry Smith   Vec              local,global;
24547c6ae99SBarry Smith   VecScatter       ltog,gtol;
246*45b6f7e9SBarry Smith   IS               to,from;
24747c6ae99SBarry Smith   PetscErrorCode   ierr;
24847c6ae99SBarry Smith 
24947c6ae99SBarry Smith   PetscFunctionBegin;
250bff4a2f0SMatthew 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");
25147c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
2523855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES)
2533855c12bSBarry 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);
2543855c12bSBarry Smith #endif
25547c6ae99SBarry Smith 
25647c6ae99SBarry Smith   if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
25747c6ae99SBarry Smith   if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
25847c6ae99SBarry Smith 
25947c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
26047c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
26147c6ae99SBarry Smith 
26247c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
26347c6ae99SBarry Smith     if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
26447c6ae99SBarry Smith     else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
26547c6ae99SBarry Smith   }
26647c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
26747c6ae99SBarry Smith     if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
26847c6ae99SBarry Smith     else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
26947c6ae99SBarry Smith   }
27047c6ae99SBarry Smith 
27147c6ae99SBarry Smith   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
27247c6ae99SBarry Smith     if (n != PETSC_DECIDE) {
27347c6ae99SBarry Smith       m = size/n;
27447c6ae99SBarry Smith     } else if (m != PETSC_DECIDE) {
27547c6ae99SBarry Smith       n = size/m;
27647c6ae99SBarry Smith     } else {
27747c6ae99SBarry Smith       /* try for squarish distribution */
278369cc0aeSBarry Smith       m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N)));
27947c6ae99SBarry Smith       if (!m) m = 1;
28047c6ae99SBarry Smith       while (m > 0) {
28147c6ae99SBarry Smith         n = size/m;
28247c6ae99SBarry Smith         if (m*n == size) break;
28347c6ae99SBarry Smith         m--;
28447c6ae99SBarry Smith       }
28547c6ae99SBarry Smith       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
28647c6ae99SBarry Smith     }
28747c6ae99SBarry 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 ");
28847c6ae99SBarry Smith   } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
28947c6ae99SBarry Smith 
29047c6ae99SBarry Smith   if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
29147c6ae99SBarry Smith   if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
29247c6ae99SBarry Smith 
29347c6ae99SBarry Smith   /*
29447c6ae99SBarry Smith      Determine locally owned region
29547c6ae99SBarry Smith      xs is the first local node number, x is the number of local nodes
29647c6ae99SBarry Smith   */
29747c6ae99SBarry Smith   if (!lx) {
298785e854fSJed Brown     ierr = PetscMalloc1(m, &dd->lx);CHKERRQ(ierr);
29947c6ae99SBarry Smith     lx   = dd->lx;
30047c6ae99SBarry Smith     for (i=0; i<m; i++) {
30147c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > i);
30247c6ae99SBarry Smith     }
30347c6ae99SBarry Smith   }
30447c6ae99SBarry Smith   x  = lx[rank % m];
30547c6ae99SBarry Smith   xs = 0;
30647c6ae99SBarry Smith   for (i=0; i<(rank % m); i++) {
30747c6ae99SBarry Smith     xs += lx[i];
30847c6ae99SBarry Smith   }
30947c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
31047c6ae99SBarry Smith   left = xs;
31147c6ae99SBarry Smith   for (i=(rank % m); i<m; i++) {
31247c6ae99SBarry Smith     left += lx[i];
31347c6ae99SBarry Smith   }
31447c6ae99SBarry 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);
31547c6ae99SBarry Smith #endif
31647c6ae99SBarry Smith 
31747c6ae99SBarry Smith   /*
31847c6ae99SBarry Smith      Determine locally owned region
31947c6ae99SBarry Smith      ys is the first local node number, y is the number of local nodes
32047c6ae99SBarry Smith   */
32147c6ae99SBarry Smith   if (!ly) {
322785e854fSJed Brown     ierr = PetscMalloc1(n, &dd->ly);CHKERRQ(ierr);
32347c6ae99SBarry Smith     ly   = dd->ly;
32447c6ae99SBarry Smith     for (i=0; i<n; i++) {
32547c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > i);
32647c6ae99SBarry Smith     }
32747c6ae99SBarry Smith   }
32847c6ae99SBarry Smith   y  = ly[rank/m];
32947c6ae99SBarry Smith   ys = 0;
33047c6ae99SBarry Smith   for (i=0; i<(rank/m); i++) {
33147c6ae99SBarry Smith     ys += ly[i];
33247c6ae99SBarry Smith   }
33347c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
33447c6ae99SBarry Smith   left = ys;
33547c6ae99SBarry Smith   for (i=(rank/m); i<n; i++) {
33647c6ae99SBarry Smith     left += ly[i];
33747c6ae99SBarry Smith   }
33847c6ae99SBarry 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);
33947c6ae99SBarry Smith #endif
34047c6ae99SBarry Smith 
341bcea557cSEthan Coon   /*
342bcea557cSEthan Coon    check if the scatter requires more than one process neighbor or wraps around
343bcea557cSEthan Coon    the domain more than once
344bcea557cSEthan Coon   */
345bff4a2f0SMatthew 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);
346bff4a2f0SMatthew 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);
34747c6ae99SBarry Smith   xe = xs + x;
34847c6ae99SBarry Smith   ye = ys + y;
34947c6ae99SBarry Smith 
350ce00eea3SSatish Balay   /* determine ghost region (Xs) and region scattered into (IXs)  */
351d9c9ebe5SPeter Brune   if (xs-s > 0) {
352d9c9ebe5SPeter Brune     Xs = xs - s; IXs = xs - s;
35388661749SPeter Brune   } else {
35488661749SPeter Brune     if (bx) {
35588661749SPeter Brune       Xs = xs - s;
35688661749SPeter Brune     } else {
35788661749SPeter Brune       Xs = 0;
35888661749SPeter Brune     }
35988661749SPeter Brune     IXs = 0;
36088661749SPeter Brune   }
361d9c9ebe5SPeter Brune   if (xe+s <= M) {
362d9c9ebe5SPeter Brune     Xe = xe + s; IXe = xe + s;
36388661749SPeter Brune   } else {
36488661749SPeter Brune     if (bx) {
365d9c9ebe5SPeter Brune       Xs = xs - s; Xe = xe + s;
36688661749SPeter Brune     } else {
36788661749SPeter Brune       Xe = M;
36888661749SPeter Brune     }
36988661749SPeter Brune     IXe = M;
37088661749SPeter Brune   }
37147c6ae99SBarry Smith 
372bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) {
373d9c9ebe5SPeter Brune     IXs = xs - s;
374d9c9ebe5SPeter Brune     IXe = xe + s;
375d9c9ebe5SPeter Brune     Xs  = xs - s;
376d9c9ebe5SPeter Brune     Xe  = xe + s;
37788661749SPeter Brune   }
37847c6ae99SBarry Smith 
379d9c9ebe5SPeter Brune   if (ys-s > 0) {
380d9c9ebe5SPeter Brune     Ys = ys - s; IYs = ys - s;
38188661749SPeter Brune   } else {
38288661749SPeter Brune     if (by) {
38388661749SPeter Brune       Ys = ys - s;
38488661749SPeter Brune     } else {
38588661749SPeter Brune       Ys = 0;
38688661749SPeter Brune     }
38788661749SPeter Brune     IYs = 0;
38888661749SPeter Brune   }
389d9c9ebe5SPeter Brune   if (ye+s <= N) {
390d9c9ebe5SPeter Brune     Ye = ye + s; IYe = ye + s;
39188661749SPeter Brune   } else {
39288661749SPeter Brune     if (by) {
39388661749SPeter Brune       Ye = ye + s;
39488661749SPeter Brune     } else {
39588661749SPeter Brune       Ye = N;
39688661749SPeter Brune     }
39788661749SPeter Brune     IYe = N;
39888661749SPeter Brune   }
39988661749SPeter Brune 
400bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) {
401d9c9ebe5SPeter Brune     IYs = ys - s;
402d9c9ebe5SPeter Brune     IYe = ye + s;
403d9c9ebe5SPeter Brune     Ys  = ys - s;
404d9c9ebe5SPeter Brune     Ye  = ye + s;
40588661749SPeter Brune   }
40688661749SPeter Brune 
40788661749SPeter Brune   /* stencil length in each direction */
408d9c9ebe5SPeter Brune   s_x = s;
409d9c9ebe5SPeter Brune   s_y = s;
41047c6ae99SBarry Smith 
41147c6ae99SBarry Smith   /* determine starting point of each processor */
41247c6ae99SBarry Smith   nn       = x*y;
413dcca6d9dSJed Brown   ierr     = PetscMalloc2(size+1,&bases,size,&ldims);CHKERRQ(ierr);
41447c6ae99SBarry Smith   ierr     = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
41547c6ae99SBarry Smith   bases[0] = 0;
41647c6ae99SBarry Smith   for (i=1; i<=size; i++) {
41747c6ae99SBarry Smith     bases[i] = ldims[i-1];
41847c6ae99SBarry Smith   }
41947c6ae99SBarry Smith   for (i=1; i<=size; i++) {
42047c6ae99SBarry Smith     bases[i] += bases[i-1];
42147c6ae99SBarry Smith   }
422ce00eea3SSatish Balay   base = bases[rank]*dof;
42347c6ae99SBarry Smith 
42447c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
425ce00eea3SSatish Balay   dd->Nlocal = x*y*dof;
426778a2246SBarry Smith   ierr       = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr);
427ce00eea3SSatish Balay   dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof;
428778a2246SBarry Smith   ierr       = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,0,&local);CHKERRQ(ierr);
42947c6ae99SBarry Smith 
43047c6ae99SBarry Smith   /* generate appropriate vector scatters */
43147c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
43247c6ae99SBarry Smith   ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr);
433ce00eea3SSatish Balay   ierr = ISCreateStride(comm,x*y*dof,start,1,&to);CHKERRQ(ierr);
43447c6ae99SBarry Smith 
435785e854fSJed Brown   ierr  = PetscMalloc1(x*y,&idx);CHKERRQ(ierr);
436ce00eea3SSatish Balay   left  = xs - Xs; right = left + x;
437ce00eea3SSatish Balay   down  = ys - Ys; up = down + y;
43847c6ae99SBarry Smith   count = 0;
43947c6ae99SBarry Smith   for (i=down; i<up; i++) {
440ce00eea3SSatish Balay     for (j=left; j<right; j++) {
441ce00eea3SSatish Balay       idx[count++] = i*(Xe-Xs) + j;
44247c6ae99SBarry Smith     }
44347c6ae99SBarry Smith   }
44447c6ae99SBarry Smith 
445ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
44647c6ae99SBarry Smith   ierr = VecScatterCreate(local,from,global,to,&ltog);CHKERRQ(ierr);
447f4a8f317SJed Brown   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)ltog);CHKERRQ(ierr);
448fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
449fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
45047c6ae99SBarry Smith 
451ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
452ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
453d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_BOX) {
454db87c5ecSEthan Coon     count = (IXe-IXs)*(IYe-IYs);
455785e854fSJed Brown     ierr  = PetscMalloc1(count,&idx);CHKERRQ(ierr);
456ce00eea3SSatish Balay 
457ce00eea3SSatish Balay     left  = IXs - Xs; right = left + (IXe-IXs);
458ce00eea3SSatish Balay     down  = IYs - Ys; up = down + (IYe-IYs);
459ce00eea3SSatish Balay     count = 0;
460ce00eea3SSatish Balay     for (i=down; i<up; i++) {
461ce00eea3SSatish Balay       for (j=left; j<right; j++) {
462ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
463ce00eea3SSatish Balay       }
464ce00eea3SSatish Balay     }
465ce00eea3SSatish Balay     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
466ce00eea3SSatish Balay 
46747c6ae99SBarry Smith   } else {
46847c6ae99SBarry Smith     /* must drop into cross shape region */
46947c6ae99SBarry Smith     /*       ---------|
47047c6ae99SBarry Smith             |  top    |
471ce00eea3SSatish Balay          |---         ---| up
47247c6ae99SBarry Smith          |   middle      |
47347c6ae99SBarry Smith          |               |
474ce00eea3SSatish Balay          ----         ---- down
47547c6ae99SBarry Smith             | bottom  |
47647c6ae99SBarry Smith             -----------
47747c6ae99SBarry Smith          Xs xs        xe Xe */
478db87c5ecSEthan Coon     count = (ys-IYs)*x + y*(IXe-IXs) + (IYe-ye)*x;
479785e854fSJed Brown     ierr  = PetscMalloc1(count,&idx);CHKERRQ(ierr);
480ce00eea3SSatish Balay 
481ce00eea3SSatish Balay     left  = xs - Xs; right = left + x;
482ce00eea3SSatish Balay     down  = ys - Ys; up = down + y;
48347c6ae99SBarry Smith     count = 0;
484ce00eea3SSatish Balay     /* bottom */
485ce00eea3SSatish Balay     for (i=(IYs-Ys); i<down; i++) {
486ce00eea3SSatish Balay       for (j=left; j<right; j++) {
487ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
48847c6ae99SBarry Smith       }
48947c6ae99SBarry Smith     }
49047c6ae99SBarry Smith     /* middle */
49147c6ae99SBarry Smith     for (i=down; i<up; i++) {
492ce00eea3SSatish Balay       for (j=(IXs-Xs); j<(IXe-Xs); j++) {
493ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
49447c6ae99SBarry Smith       }
49547c6ae99SBarry Smith     }
49647c6ae99SBarry Smith     /* top */
497ce00eea3SSatish Balay     for (i=up; i<up+IYe-ye; i++) {
498ce00eea3SSatish Balay       for (j=left; j<right; j++) {
499ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
50047c6ae99SBarry Smith       }
50147c6ae99SBarry Smith     }
50247c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
50347c6ae99SBarry Smith   }
50447c6ae99SBarry Smith 
50547c6ae99SBarry Smith 
50647c6ae99SBarry Smith   /* determine who lies on each side of us stored in    n6 n7 n8
50747c6ae99SBarry Smith                                                         n3    n5
50847c6ae99SBarry Smith                                                         n0 n1 n2
50947c6ae99SBarry Smith   */
51047c6ae99SBarry Smith 
51147c6ae99SBarry Smith   /* Assume the Non-Periodic Case */
51247c6ae99SBarry Smith   n1 = rank - m;
51347c6ae99SBarry Smith   if (rank % m) {
51447c6ae99SBarry Smith     n0 = n1 - 1;
51547c6ae99SBarry Smith   } else {
51647c6ae99SBarry Smith     n0 = -1;
51747c6ae99SBarry Smith   }
51847c6ae99SBarry Smith   if ((rank+1) % m) {
51947c6ae99SBarry Smith     n2 = n1 + 1;
52047c6ae99SBarry Smith     n5 = rank + 1;
52147c6ae99SBarry Smith     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
52247c6ae99SBarry Smith   } else {
52347c6ae99SBarry Smith     n2 = -1; n5 = -1; n8 = -1;
52447c6ae99SBarry Smith   }
52547c6ae99SBarry Smith   if (rank % m) {
52647c6ae99SBarry Smith     n3 = rank - 1;
52747c6ae99SBarry Smith     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
52847c6ae99SBarry Smith   } else {
52947c6ae99SBarry Smith     n3 = -1; n6 = -1;
53047c6ae99SBarry Smith   }
53147c6ae99SBarry Smith   n7 = rank + m; if (n7 >= m*n) n7 = -1;
53247c6ae99SBarry Smith 
533bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC && by == DM_BOUNDARY_PERIODIC) {
53447c6ae99SBarry Smith     /* Modify for Periodic Cases */
53547c6ae99SBarry Smith     /* Handle all four corners */
53647c6ae99SBarry Smith     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
53747c6ae99SBarry Smith     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
53847c6ae99SBarry Smith     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
53947c6ae99SBarry Smith     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
54047c6ae99SBarry Smith 
54147c6ae99SBarry Smith     /* Handle Top and Bottom Sides */
54247c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
54347c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
54447c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
54547c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
54647c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
54747c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
54847c6ae99SBarry Smith 
54947c6ae99SBarry Smith     /* Handle Left and Right Sides */
55047c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
55147c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
55247c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
55347c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
55447c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
55547c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
556bff4a2f0SMatthew G. Knepley   } else if (by == DM_BOUNDARY_PERIODIC) {  /* Handle Top and Bottom Sides */
557ce00eea3SSatish Balay     if (n1 < 0) n1 = rank + m * (n-1);
558ce00eea3SSatish Balay     if (n7 < 0) n7 = rank - m * (n-1);
559ce00eea3SSatish Balay     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
560ce00eea3SSatish Balay     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
561ce00eea3SSatish Balay     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
562ce00eea3SSatish Balay     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
563bff4a2f0SMatthew G. Knepley   } else if (bx == DM_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */
564ce00eea3SSatish Balay     if (n3 < 0) n3 = rank + (m-1);
565ce00eea3SSatish Balay     if (n5 < 0) n5 = rank - (m-1);
566ce00eea3SSatish Balay     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
567ce00eea3SSatish Balay     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
568ce00eea3SSatish Balay     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
569ce00eea3SSatish Balay     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
57047c6ae99SBarry Smith   }
571ce00eea3SSatish Balay 
572785e854fSJed Brown   ierr = PetscMalloc1(9,&dd->neighbors);CHKERRQ(ierr);
5738865f1eaSKarl Rupp 
57447c6ae99SBarry Smith   dd->neighbors[0] = n0;
57547c6ae99SBarry Smith   dd->neighbors[1] = n1;
57647c6ae99SBarry Smith   dd->neighbors[2] = n2;
57747c6ae99SBarry Smith   dd->neighbors[3] = n3;
57847c6ae99SBarry Smith   dd->neighbors[4] = rank;
57947c6ae99SBarry Smith   dd->neighbors[5] = n5;
58047c6ae99SBarry Smith   dd->neighbors[6] = n6;
58147c6ae99SBarry Smith   dd->neighbors[7] = n7;
58247c6ae99SBarry Smith   dd->neighbors[8] = n8;
58347c6ae99SBarry Smith 
584d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
58547c6ae99SBarry Smith     /* save corner processor numbers */
58647c6ae99SBarry Smith     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
58747c6ae99SBarry Smith     n0  = n2 = n6 = n8 = -1;
58847c6ae99SBarry Smith   }
58947c6ae99SBarry Smith 
590785e854fSJed Brown   ierr = PetscMalloc1((Xe-Xs)*(Ye-Ys),&idx);CHKERRQ(ierr);
5913bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)da,(Xe-Xs)*(Ye-Ys)*sizeof(PetscInt));CHKERRQ(ierr);
59247c6ae99SBarry Smith 
593ce00eea3SSatish Balay   nn = 0;
59447c6ae99SBarry Smith   xbase = bases[rank];
59547c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
59647c6ae99SBarry Smith     if (n0 >= 0) { /* left below */
597ce00eea3SSatish Balay       x_t = lx[n0 % m];
59847c6ae99SBarry Smith       y_t = ly[(n0/m)];
59947c6ae99SBarry Smith       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
6008865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
60147c6ae99SBarry Smith     }
602ac119b13SBarry Smith 
60347c6ae99SBarry Smith     if (n1 >= 0) { /* directly below */
60447c6ae99SBarry Smith       x_t = x;
60547c6ae99SBarry Smith       y_t = ly[(n1/m)];
60647c6ae99SBarry Smith       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
6078865f1eaSKarl Rupp       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
608bff4a2f0SMatthew G. Knepley     } else if (by == DM_BOUNDARY_MIRROR) {
6098865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
61047c6ae99SBarry Smith     }
611ac119b13SBarry Smith 
61247c6ae99SBarry Smith     if (n2 >= 0) { /* right below */
613ce00eea3SSatish Balay       x_t = lx[n2 % m];
61447c6ae99SBarry Smith       y_t = ly[(n2/m)];
61547c6ae99SBarry Smith       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
6168865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
61747c6ae99SBarry Smith     }
61847c6ae99SBarry Smith   }
61947c6ae99SBarry Smith 
62047c6ae99SBarry Smith   for (i=0; i<y; i++) {
62147c6ae99SBarry Smith     if (n3 >= 0) { /* directly left */
622ce00eea3SSatish Balay       x_t = lx[n3 % m];
62347c6ae99SBarry Smith       /* y_t = y; */
62447c6ae99SBarry Smith       s_t = bases[n3] + (i+1)*x_t - s_x;
6258865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
626bff4a2f0SMatthew G. Knepley     } else if (bx == DM_BOUNDARY_MIRROR) {
6278865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
62847c6ae99SBarry Smith     }
62947c6ae99SBarry Smith 
6308865f1eaSKarl Rupp     for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
63147c6ae99SBarry Smith 
63247c6ae99SBarry Smith     if (n5 >= 0) { /* directly right */
633ce00eea3SSatish Balay       x_t = lx[n5 % m];
63447c6ae99SBarry Smith       /* y_t = y; */
63547c6ae99SBarry Smith       s_t = bases[n5] + (i)*x_t;
6368865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
637bff4a2f0SMatthew G. Knepley     } else if (bx == DM_BOUNDARY_MIRROR) {
6388865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
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;
6478865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
64847c6ae99SBarry Smith     }
649ac119b13SBarry Smith 
65047c6ae99SBarry Smith     if (n7 >= 0) { /* directly above */
65147c6ae99SBarry Smith       x_t = x;
65247c6ae99SBarry Smith       /* y_t = ly[(n7/m)]; */
65347c6ae99SBarry Smith       s_t = bases[n7] + (i-1)*x_t;
6548865f1eaSKarl Rupp       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
655bff4a2f0SMatthew G. Knepley     } else if (by == DM_BOUNDARY_MIRROR) {
6568865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
65747c6ae99SBarry Smith     }
658ac119b13SBarry Smith 
65947c6ae99SBarry Smith     if (n8 >= 0) { /* right above */
660ce00eea3SSatish Balay       x_t = lx[n8 % m];
66147c6ae99SBarry Smith       /* y_t = ly[(n8/m)]; */
66247c6ae99SBarry Smith       s_t = bases[n8] + (i-1)*x_t;
6638865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
66447c6ae99SBarry Smith     }
66547c6ae99SBarry Smith   }
66647c6ae99SBarry Smith 
667ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
66847c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
6693bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr);
670fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
671fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
67247c6ae99SBarry Smith 
673d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
674ce00eea3SSatish Balay     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
675ce00eea3SSatish Balay   }
676ce00eea3SSatish Balay 
67788661749SPeter Brune   if (((stencil_type == DMDA_STENCIL_STAR)  ||
678bff4a2f0SMatthew G. Knepley        (bx && bx != DM_BOUNDARY_PERIODIC) ||
679bff4a2f0SMatthew G. Knepley        (by && by != DM_BOUNDARY_PERIODIC))) {
68047c6ae99SBarry Smith     /*
68147c6ae99SBarry Smith         Recompute the local to global mappings, this time keeping the
682ce00eea3SSatish Balay       information about the cross corner processor numbers and any ghosted
683ce00eea3SSatish Balay       but not periodic indices.
68447c6ae99SBarry Smith     */
68547c6ae99SBarry Smith     nn    = 0;
68647c6ae99SBarry Smith     xbase = bases[rank];
68747c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
68847c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
689ce00eea3SSatish Balay         x_t = lx[n0 % m];
69047c6ae99SBarry Smith         y_t = ly[(n0/m)];
69147c6ae99SBarry Smith         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
6928865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
693ce00eea3SSatish Balay       } else if (xs-Xs > 0 && ys-Ys > 0) {
6948865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
69547c6ae99SBarry Smith       }
69647c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
69747c6ae99SBarry Smith         x_t = x;
69847c6ae99SBarry Smith         y_t = ly[(n1/m)];
69947c6ae99SBarry Smith         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
7008865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
701ce00eea3SSatish Balay       } else if (ys-Ys > 0) {
702bff4a2f0SMatthew G. Knepley         if (by == DM_BOUNDARY_MIRROR) {
7038865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
704624904c4SBarry Smith         } else {
7058865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
70647c6ae99SBarry Smith         }
707624904c4SBarry Smith       }
70847c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
709ce00eea3SSatish Balay         x_t = lx[n2 % m];
71047c6ae99SBarry Smith         y_t = ly[(n2/m)];
71147c6ae99SBarry Smith         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
7128865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
713ce00eea3SSatish Balay       } else if (Xe-xe> 0 && ys-Ys > 0) {
7148865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
71547c6ae99SBarry Smith       }
71647c6ae99SBarry Smith     }
71747c6ae99SBarry Smith 
71847c6ae99SBarry Smith     for (i=0; i<y; i++) {
71947c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
720ce00eea3SSatish Balay         x_t = lx[n3 % m];
72147c6ae99SBarry Smith         /* y_t = y; */
72247c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x;
7238865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
724ce00eea3SSatish Balay       } else if (xs-Xs > 0) {
725bff4a2f0SMatthew G. Knepley         if (bx == DM_BOUNDARY_MIRROR) {
7268865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
727624904c4SBarry Smith         } else {
7288865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
72947c6ae99SBarry Smith         }
730624904c4SBarry Smith       }
73147c6ae99SBarry Smith 
7328865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
73347c6ae99SBarry Smith 
73447c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
735ce00eea3SSatish Balay         x_t = lx[n5 % m];
73647c6ae99SBarry Smith         /* y_t = y; */
73747c6ae99SBarry Smith         s_t = bases[n5] + (i)*x_t;
7388865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
739ce00eea3SSatish Balay       } else if (Xe-xe > 0) {
740bff4a2f0SMatthew G. Knepley         if (bx == DM_BOUNDARY_MIRROR) {
7418865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
742624904c4SBarry Smith         } else {
7438865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
74447c6ae99SBarry Smith         }
74547c6ae99SBarry Smith       }
746624904c4SBarry Smith     }
74747c6ae99SBarry Smith 
74847c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
74947c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
750ce00eea3SSatish Balay         x_t = lx[n6 % m];
75147c6ae99SBarry Smith         /* y_t = ly[(n6/m)]; */
75247c6ae99SBarry Smith         s_t = bases[n6] + (i)*x_t - s_x;
7538865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
754ce00eea3SSatish Balay       } else if (xs-Xs > 0 && Ye-ye > 0) {
7558865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
75647c6ae99SBarry Smith       }
75747c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
75847c6ae99SBarry Smith         x_t = x;
75947c6ae99SBarry Smith         /* y_t = ly[(n7/m)]; */
76047c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t;
7618865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
762ce00eea3SSatish Balay       } else if (Ye-ye > 0) {
763bff4a2f0SMatthew G. Knepley         if (by == DM_BOUNDARY_MIRROR) {
7648865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
765624904c4SBarry Smith         } else {
7668865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
76747c6ae99SBarry Smith         }
768624904c4SBarry Smith       }
76947c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
770ce00eea3SSatish Balay         x_t = lx[n8 % m];
77147c6ae99SBarry Smith         /* y_t = ly[(n8/m)]; */
77247c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t;
7738865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
774ce00eea3SSatish Balay       } else if (Xe-xe > 0 && Ye-ye > 0) {
7758865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
77647c6ae99SBarry Smith       }
77747c6ae99SBarry Smith     }
77847c6ae99SBarry Smith   }
779ce00eea3SSatish Balay   /*
780ce00eea3SSatish Balay      Set the local to global ordering in the global vector, this allows use
781ce00eea3SSatish Balay      of VecSetValuesLocal().
782ce00eea3SSatish Balay   */
783*45b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr);
7843bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr);
78547c6ae99SBarry Smith 
786ce00eea3SSatish Balay   ierr  = PetscFree2(bases,ldims);CHKERRQ(ierr);
78747c6ae99SBarry Smith   dd->m = m;  dd->n  = n;
788ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
789ce00eea3SSatish Balay   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
790ce00eea3SSatish Balay   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
79147c6ae99SBarry Smith 
792fcfd50ebSBarry Smith   ierr = VecDestroy(&local);CHKERRQ(ierr);
793fcfd50ebSBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
79447c6ae99SBarry Smith 
79547c6ae99SBarry Smith   dd->gtol      = gtol;
79647c6ae99SBarry Smith   dd->ltog      = ltog;
79747c6ae99SBarry Smith   dd->base      = base;
7989a42bb27SBarry Smith   da->ops->view = DMView_DA_2d;
7990298fd71SBarry Smith   dd->ltol      = NULL;
8000298fd71SBarry Smith   dd->ao        = NULL;
80147c6ae99SBarry Smith   PetscFunctionReturn(0);
80247c6ae99SBarry Smith }
80347c6ae99SBarry Smith 
80447c6ae99SBarry Smith #undef __FUNCT__
805aa219208SBarry Smith #define __FUNCT__ "DMDACreate2d"
80647c6ae99SBarry Smith /*@C
807aa219208SBarry Smith    DMDACreate2d -  Creates an object that will manage the communication of  two-dimensional
80847c6ae99SBarry Smith    regular array data that is distributed across some processors.
80947c6ae99SBarry Smith 
81047c6ae99SBarry Smith    Collective on MPI_Comm
81147c6ae99SBarry Smith 
81247c6ae99SBarry Smith    Input Parameters:
81347c6ae99SBarry Smith +  comm - MPI communicator
8141321219cSEthan Coon .  bx,by - type of ghost nodes the array have.
815bff4a2f0SMatthew G. Knepley          Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC.
816aa219208SBarry Smith .  stencil_type - stencil type.  Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
81747c6ae99SBarry 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
81847c6ae99SBarry Smith             from the command line with -da_grid_x <M> -da_grid_y <N>)
81947c6ae99SBarry Smith .  m,n - corresponding number of processors in each dimension
82047c6ae99SBarry Smith          (or PETSC_DECIDE to have calculated)
82147c6ae99SBarry Smith .  dof - number of degrees of freedom per node
82247c6ae99SBarry Smith .  s - stencil width
82347c6ae99SBarry Smith -  lx, ly - arrays containing the number of nodes in each cell along
8240298fd71SBarry Smith            the x and y coordinates, or NULL. If non-null, these
82547c6ae99SBarry Smith            must be of length as m and n, and the corresponding
82647c6ae99SBarry Smith            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
82747c6ae99SBarry Smith            must be M, and the sum of the ly[] entries must be N.
82847c6ae99SBarry Smith 
82947c6ae99SBarry Smith    Output Parameter:
83047c6ae99SBarry Smith .  da - the resulting distributed array object
83147c6ae99SBarry Smith 
83247c6ae99SBarry Smith    Options Database Key:
833706a11cbSBarry Smith +  -dm_view - Calls DMView() at the conclusion of DMDACreate2d()
83447c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
83547c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
83647c6ae99SBarry Smith .  -da_processors_x <nx> - number of processors in x direction
83747c6ae99SBarry Smith .  -da_processors_y <ny> - number of processors in y direction
838e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement ratio in x direction
839e0f5d30fSBarry Smith .  -da_refine_y <ry> - refinement ratio in y direction
840e0f5d30fSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating, if M or N < 0
841e0f5d30fSBarry Smith 
84247c6ae99SBarry Smith 
84347c6ae99SBarry Smith    Level: beginner
84447c6ae99SBarry Smith 
84547c6ae99SBarry Smith    Notes:
846aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
847aa219208SBarry Smith    standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
84847c6ae99SBarry Smith    the standard 9-pt stencil.
84947c6ae99SBarry Smith 
850aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
851564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
852564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
85347c6ae99SBarry Smith 
85447c6ae99SBarry Smith .keywords: distributed array, create, two-dimensional
85547c6ae99SBarry Smith 
856aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
85799f0b487SRichard Tran Mills           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(),
858d461ba97SBarry Smith           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
85947c6ae99SBarry Smith 
86047c6ae99SBarry Smith @*/
861fe16a2e9SBarry Smith 
862bff4a2f0SMatthew G. Knepley PetscErrorCode  DMDACreate2d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMDAStencilType stencil_type,
8639a42bb27SBarry Smith                           PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
86447c6ae99SBarry Smith {
86547c6ae99SBarry Smith   PetscErrorCode ierr;
86647c6ae99SBarry Smith 
86747c6ae99SBarry Smith   PetscFunctionBegin;
868aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
869aa219208SBarry Smith   ierr = DMDASetDim(*da, 2);CHKERRQ(ierr);
870aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, 1);CHKERRQ(ierr);
871aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr);
872bff4a2f0SMatthew G. Knepley   ierr = DMDASetBoundaryType(*da, bx, by, DM_BOUNDARY_NONE);CHKERRQ(ierr);
873aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
874aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
875aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
8760298fd71SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, NULL);CHKERRQ(ierr);
87747c6ae99SBarry Smith   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
8789a42bb27SBarry Smith   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
8799a42bb27SBarry Smith   ierr = DMSetUp(*da);CHKERRQ(ierr);
88047c6ae99SBarry Smith   PetscFunctionReturn(0);
88147c6ae99SBarry Smith }
882