xref: /petsc/src/dm/impls/da/da2.c (revision 5b399a637f946da949d4901b9cfb4dd7404a5662)
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);
52*5b399a63SLisandro Dalcin     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
53*5b399a63SLisandro Dalcin     if (isnull) PetscFunctionReturn(0);
5447c6ae99SBarry Smith 
55*5b399a63SLisandro Dalcin     ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
56*5b399a63SLisandro Dalcin     ierr = PetscDrawClear(draw);CHKERRQ(ierr);
57*5b399a63SLisandro Dalcin 
58*5b399a63SLisandro Dalcin     ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
59*5b399a63SLisandro Dalcin 
60*5b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
6147c6ae99SBarry Smith     /* first processor draw all node lines */
6247c6ae99SBarry Smith     if (!rank) {
6347c6ae99SBarry Smith       ymin = 0.0; ymax = dd->N - 1;
6447c6ae99SBarry Smith       for (xmin=0; xmin<dd->M; xmin++) {
6547c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
6647c6ae99SBarry Smith       }
6747c6ae99SBarry Smith       xmin = 0.0; xmax = dd->M - 1;
6847c6ae99SBarry Smith       for (ymin=0; ymin<dd->N; ymin++) {
6947c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
7047c6ae99SBarry Smith       }
7147c6ae99SBarry Smith     }
72*5b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
73*5b399a63SLisandro Dalcin 
74*5b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
7547c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
7647c6ae99SBarry Smith 
77*5b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
7847c6ae99SBarry Smith     /* draw my box */
79*5b399a63SLisandro Dalcin     xmin = dd->xs/dd->w; xmax =(dd->xe-1)/dd->w; ymin = dd->ys; ymax = dd->ye - 1;
8047c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
8147c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
8247c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
8347c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
8447c6ae99SBarry Smith     /* put in numbers */
8547c6ae99SBarry Smith     base = (dd->base)/dd->w;
8647c6ae99SBarry Smith     for (y=ymin; y<=ymax; y++) {
8747c6ae99SBarry Smith       for (x=xmin; x<=xmax; x++) {
88*5b399a63SLisandro Dalcin         ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)base++);CHKERRQ(ierr);
8947c6ae99SBarry Smith         ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
9047c6ae99SBarry Smith       }
9147c6ae99SBarry Smith     }
92*5b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
9347c6ae99SBarry Smith 
94*5b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
9547c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
9647c6ae99SBarry Smith 
97*5b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
98*5b399a63SLisandro Dalcin     /* overlay ghost numbers, useful for error checking */
9945b6f7e9SBarry Smith     ierr = ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr);
100*5b399a63SLisandro Dalcin     base = 0; xmin = dd->Xs; xmax = dd->Xe; ymin = dd->Ys; ymax = dd->Ye;
10147c6ae99SBarry Smith     for (y=ymin; y<ymax; y++) {
10247c6ae99SBarry Smith       for (x=xmin; x<xmax; x++) {
10347c6ae99SBarry Smith         if ((base % dd->w) == 0) {
104*5b399a63SLisandro Dalcin           ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)(idx[base/dd->w]));CHKERRQ(ierr);
10547c6ae99SBarry Smith           ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
10647c6ae99SBarry Smith         }
10747c6ae99SBarry Smith         base++;
10847c6ae99SBarry Smith       }
10947c6ae99SBarry Smith     }
110302440fdSBarry Smith     ierr = ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr);
111*5b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
112*5b399a63SLisandro Dalcin 
113*5b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
11447c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
1159a42bb27SBarry Smith   } else if (isbinary) {
1169a42bb27SBarry Smith     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
1179a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1189a42bb27SBarry Smith   } else if (ismatlab) {
1199a42bb27SBarry Smith     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
1209a42bb27SBarry Smith #endif
12111aeaf0aSBarry Smith   }
12247c6ae99SBarry Smith   PetscFunctionReturn(0);
12347c6ae99SBarry Smith }
12447c6ae99SBarry Smith 
12547c6ae99SBarry Smith /*
12647c6ae99SBarry Smith       M is number of grid points
12747c6ae99SBarry Smith       m is number of processors
12847c6ae99SBarry Smith 
12947c6ae99SBarry Smith */
13047c6ae99SBarry Smith #undef __FUNCT__
131aa219208SBarry Smith #define __FUNCT__ "DMDASplitComm2d"
1327087cfbeSBarry Smith PetscErrorCode  DMDASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm)
13347c6ae99SBarry Smith {
13447c6ae99SBarry Smith   PetscErrorCode ierr;
13547c6ae99SBarry Smith   PetscInt       m,n = 0,x = 0,y = 0;
13647c6ae99SBarry Smith   PetscMPIInt    size,csize,rank;
13747c6ae99SBarry Smith 
13847c6ae99SBarry Smith   PetscFunctionBegin;
13947c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
14047c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
14147c6ae99SBarry Smith 
14247c6ae99SBarry Smith   csize = 4*size;
14347c6ae99SBarry Smith   do {
14447c6ae99SBarry 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);
14547c6ae99SBarry Smith     csize = csize/4;
14647c6ae99SBarry Smith 
147369cc0aeSBarry Smith     m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)csize)/((PetscReal)N)));
14847c6ae99SBarry Smith     if (!m) m = 1;
14947c6ae99SBarry Smith     while (m > 0) {
15047c6ae99SBarry Smith       n = csize/m;
15147c6ae99SBarry Smith       if (m*n == csize) break;
15247c6ae99SBarry Smith       m--;
15347c6ae99SBarry Smith     }
15447c6ae99SBarry Smith     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
15547c6ae99SBarry Smith 
15647c6ae99SBarry Smith     x = M/m + ((M % m) > ((csize-1) % m));
15747c6ae99SBarry Smith     y = (N + (csize-1)/m)/n;
15847c6ae99SBarry Smith   } while ((x < 4 || y < 4) && csize > 1);
15947c6ae99SBarry Smith   if (size != csize) {
16047c6ae99SBarry Smith     MPI_Group   entire_group,sub_group;
16147c6ae99SBarry Smith     PetscMPIInt i,*groupies;
16247c6ae99SBarry Smith 
16347c6ae99SBarry Smith     ierr = MPI_Comm_group(comm,&entire_group);CHKERRQ(ierr);
164785e854fSJed Brown     ierr = PetscMalloc1(csize,&groupies);CHKERRQ(ierr);
16547c6ae99SBarry Smith     for (i=0; i<csize; i++) {
16647c6ae99SBarry Smith       groupies[i] = (rank/csize)*csize + i;
16747c6ae99SBarry Smith     }
16847c6ae99SBarry Smith     ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);CHKERRQ(ierr);
16947c6ae99SBarry Smith     ierr = PetscFree(groupies);CHKERRQ(ierr);
17047c6ae99SBarry Smith     ierr = MPI_Comm_create(comm,sub_group,outcomm);CHKERRQ(ierr);
17147c6ae99SBarry Smith     ierr = MPI_Group_free(&entire_group);CHKERRQ(ierr);
17247c6ae99SBarry Smith     ierr = MPI_Group_free(&sub_group);CHKERRQ(ierr);
173aa219208SBarry Smith     ierr = PetscInfo1(0,"DMDASplitComm2d:Creating redundant coarse problems of size %d\n",csize);CHKERRQ(ierr);
17447c6ae99SBarry Smith   } else {
17547c6ae99SBarry Smith     *outcomm = comm;
17647c6ae99SBarry Smith   }
17747c6ae99SBarry Smith   PetscFunctionReturn(0);
17847c6ae99SBarry Smith }
17947c6ae99SBarry Smith 
18047c6ae99SBarry Smith #if defined(new)
18147c6ae99SBarry Smith #undef __FUNCT__
182aa219208SBarry Smith #define __FUNCT__ "DMDAGetDiagonal_MFFD"
18347c6ae99SBarry Smith /*
184aa219208SBarry Smith   DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
185aa219208SBarry Smith     function lives on a DMDA
18647c6ae99SBarry Smith 
18747c6ae99SBarry Smith         y ~= (F(u + ha) - F(u))/h,
18847c6ae99SBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
18947c6ae99SBarry Smith         u = current iterate
19047c6ae99SBarry Smith         h = difference interval
19147c6ae99SBarry Smith */
192aa219208SBarry Smith PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a)
19347c6ae99SBarry Smith {
19447c6ae99SBarry Smith   PetscScalar    h,*aa,*ww,v;
19547c6ae99SBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
19647c6ae99SBarry Smith   PetscErrorCode ierr;
19747c6ae99SBarry Smith   PetscInt       gI,nI;
19847c6ae99SBarry Smith   MatStencil     stencil;
199aa219208SBarry Smith   DMDALocalInfo  info;
20047c6ae99SBarry Smith 
20147c6ae99SBarry Smith   PetscFunctionBegin;
20247c6ae99SBarry Smith   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
20347c6ae99SBarry Smith   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
20447c6ae99SBarry Smith 
20547c6ae99SBarry Smith   ierr = VecGetArray(U,&ww);CHKERRQ(ierr);
20647c6ae99SBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
20747c6ae99SBarry Smith 
20847c6ae99SBarry Smith   nI = 0;
20947c6ae99SBarry Smith   h  = ww[gI];
21047c6ae99SBarry Smith   if (h == 0.0) h = 1.0;
21147c6ae99SBarry Smith   if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin;
21247c6ae99SBarry Smith   else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
21347c6ae99SBarry Smith   h *= epsilon;
21447c6ae99SBarry Smith 
21547c6ae99SBarry Smith   ww[gI] += h;
21647c6ae99SBarry Smith   ierr    = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
21747c6ae99SBarry Smith   aa[nI]  = (v - aa[nI])/h;
21847c6ae99SBarry Smith   ww[gI] -= h;
21947c6ae99SBarry Smith   nI++;
2208865f1eaSKarl Rupp 
22147c6ae99SBarry Smith   ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr);
22247c6ae99SBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
22347c6ae99SBarry Smith   PetscFunctionReturn(0);
22447c6ae99SBarry Smith }
22547c6ae99SBarry Smith #endif
22647c6ae99SBarry Smith 
22747c6ae99SBarry Smith #undef __FUNCT__
2289a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_2D"
2297087cfbeSBarry Smith PetscErrorCode  DMSetUp_DA_2D(DM da)
23047c6ae99SBarry Smith {
23147c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
23247c6ae99SBarry Smith   const PetscInt   M            = dd->M;
23347c6ae99SBarry Smith   const PetscInt   N            = dd->N;
23447c6ae99SBarry Smith   PetscInt         m            = dd->m;
23547c6ae99SBarry Smith   PetscInt         n            = dd->n;
23647c6ae99SBarry Smith   const PetscInt   dof          = dd->w;
23747c6ae99SBarry Smith   const PetscInt   s            = dd->s;
238bff4a2f0SMatthew G. Knepley   DMBoundaryType   bx           = dd->bx;
239bff4a2f0SMatthew G. Knepley   DMBoundaryType   by           = dd->by;
24019fd82e9SBarry Smith   DMDAStencilType  stencil_type = dd->stencil_type;
24147c6ae99SBarry Smith   PetscInt         *lx          = dd->lx;
24247c6ae99SBarry Smith   PetscInt         *ly          = dd->ly;
24347c6ae99SBarry Smith   MPI_Comm         comm;
24447c6ae99SBarry Smith   PetscMPIInt      rank,size;
245bd1fc5aeSBarry Smith   PetscInt         xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,IXs,IXe,IYs,IYe;
2468ea3bf28SBarry Smith   PetscInt         up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn;
247db87c5ecSEthan Coon   PetscInt         xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
24847c6ae99SBarry Smith   PetscInt         s_x,s_y; /* s proportionalized to w */
24947c6ae99SBarry Smith   PetscInt         sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
25047c6ae99SBarry Smith   Vec              local,global;
251bd1fc5aeSBarry Smith   VecScatter       gtol;
25245b6f7e9SBarry Smith   IS               to,from;
25347c6ae99SBarry Smith   PetscErrorCode   ierr;
25447c6ae99SBarry Smith 
25547c6ae99SBarry Smith   PetscFunctionBegin;
256bff4a2f0SMatthew 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");
25747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
2583855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES)
2593855c12bSBarry 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);
2603855c12bSBarry Smith #endif
26147c6ae99SBarry Smith 
26247c6ae99SBarry Smith   if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
26347c6ae99SBarry Smith   if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
26447c6ae99SBarry Smith 
26547c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
26647c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
26747c6ae99SBarry Smith 
2687d310018SBarry Smith   dd->p = 1;
26947c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
27047c6ae99SBarry Smith     if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
27147c6ae99SBarry Smith     else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
27247c6ae99SBarry Smith   }
27347c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
27447c6ae99SBarry Smith     if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
27547c6ae99SBarry Smith     else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
27647c6ae99SBarry Smith   }
27747c6ae99SBarry Smith 
27847c6ae99SBarry Smith   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
27947c6ae99SBarry Smith     if (n != PETSC_DECIDE) {
28047c6ae99SBarry Smith       m = size/n;
28147c6ae99SBarry Smith     } else if (m != PETSC_DECIDE) {
28247c6ae99SBarry Smith       n = size/m;
28347c6ae99SBarry Smith     } else {
28447c6ae99SBarry Smith       /* try for squarish distribution */
285369cc0aeSBarry Smith       m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N)));
28647c6ae99SBarry Smith       if (!m) m = 1;
28747c6ae99SBarry Smith       while (m > 0) {
28847c6ae99SBarry Smith         n = size/m;
28947c6ae99SBarry Smith         if (m*n == size) break;
29047c6ae99SBarry Smith         m--;
29147c6ae99SBarry Smith       }
29247c6ae99SBarry Smith       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
29347c6ae99SBarry Smith     }
29447c6ae99SBarry 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 ");
29547c6ae99SBarry Smith   } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
29647c6ae99SBarry Smith 
29747c6ae99SBarry Smith   if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
29847c6ae99SBarry Smith   if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
29947c6ae99SBarry Smith 
30047c6ae99SBarry Smith   /*
30147c6ae99SBarry Smith      Determine locally owned region
30247c6ae99SBarry Smith      xs is the first local node number, x is the number of local nodes
30347c6ae99SBarry Smith   */
30447c6ae99SBarry Smith   if (!lx) {
305785e854fSJed Brown     ierr = PetscMalloc1(m, &dd->lx);CHKERRQ(ierr);
30647c6ae99SBarry Smith     lx   = dd->lx;
30747c6ae99SBarry Smith     for (i=0; i<m; i++) {
30847c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > i);
30947c6ae99SBarry Smith     }
31047c6ae99SBarry Smith   }
31147c6ae99SBarry Smith   x  = lx[rank % m];
31247c6ae99SBarry Smith   xs = 0;
31347c6ae99SBarry Smith   for (i=0; i<(rank % m); i++) {
31447c6ae99SBarry Smith     xs += lx[i];
31547c6ae99SBarry Smith   }
31647c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
31747c6ae99SBarry Smith   left = xs;
31847c6ae99SBarry Smith   for (i=(rank % m); i<m; i++) {
31947c6ae99SBarry Smith     left += lx[i];
32047c6ae99SBarry Smith   }
32147c6ae99SBarry 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);
32247c6ae99SBarry Smith #endif
32347c6ae99SBarry Smith 
32447c6ae99SBarry Smith   /*
32547c6ae99SBarry Smith      Determine locally owned region
32647c6ae99SBarry Smith      ys is the first local node number, y is the number of local nodes
32747c6ae99SBarry Smith   */
32847c6ae99SBarry Smith   if (!ly) {
329785e854fSJed Brown     ierr = PetscMalloc1(n, &dd->ly);CHKERRQ(ierr);
33047c6ae99SBarry Smith     ly   = dd->ly;
33147c6ae99SBarry Smith     for (i=0; i<n; i++) {
33247c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > i);
33347c6ae99SBarry Smith     }
33447c6ae99SBarry Smith   }
33547c6ae99SBarry Smith   y  = ly[rank/m];
33647c6ae99SBarry Smith   ys = 0;
33747c6ae99SBarry Smith   for (i=0; i<(rank/m); i++) {
33847c6ae99SBarry Smith     ys += ly[i];
33947c6ae99SBarry Smith   }
34047c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
34147c6ae99SBarry Smith   left = ys;
34247c6ae99SBarry Smith   for (i=(rank/m); i<n; i++) {
34347c6ae99SBarry Smith     left += ly[i];
34447c6ae99SBarry Smith   }
34547c6ae99SBarry 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);
34647c6ae99SBarry Smith #endif
34747c6ae99SBarry Smith 
348bcea557cSEthan Coon   /*
349bcea557cSEthan Coon    check if the scatter requires more than one process neighbor or wraps around
350bcea557cSEthan Coon    the domain more than once
351bcea557cSEthan Coon   */
352bff4a2f0SMatthew 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);
353bff4a2f0SMatthew 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);
35447c6ae99SBarry Smith   xe = xs + x;
35547c6ae99SBarry Smith   ye = ys + y;
35647c6ae99SBarry Smith 
357ce00eea3SSatish Balay   /* determine ghost region (Xs) and region scattered into (IXs)  */
358d9c9ebe5SPeter Brune   if (xs-s > 0) {
359d9c9ebe5SPeter Brune     Xs = xs - s; IXs = xs - s;
36088661749SPeter Brune   } else {
36188661749SPeter Brune     if (bx) {
36288661749SPeter Brune       Xs = xs - s;
36388661749SPeter Brune     } else {
36488661749SPeter Brune       Xs = 0;
36588661749SPeter Brune     }
36688661749SPeter Brune     IXs = 0;
36788661749SPeter Brune   }
368d9c9ebe5SPeter Brune   if (xe+s <= M) {
369d9c9ebe5SPeter Brune     Xe = xe + s; IXe = xe + s;
37088661749SPeter Brune   } else {
37188661749SPeter Brune     if (bx) {
372d9c9ebe5SPeter Brune       Xs = xs - s; Xe = xe + s;
37388661749SPeter Brune     } else {
37488661749SPeter Brune       Xe = M;
37588661749SPeter Brune     }
37688661749SPeter Brune     IXe = M;
37788661749SPeter Brune   }
37847c6ae99SBarry Smith 
379bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) {
380d9c9ebe5SPeter Brune     IXs = xs - s;
381d9c9ebe5SPeter Brune     IXe = xe + s;
382d9c9ebe5SPeter Brune     Xs  = xs - s;
383d9c9ebe5SPeter Brune     Xe  = xe + s;
38488661749SPeter Brune   }
38547c6ae99SBarry Smith 
386d9c9ebe5SPeter Brune   if (ys-s > 0) {
387d9c9ebe5SPeter Brune     Ys = ys - s; IYs = ys - s;
38888661749SPeter Brune   } else {
38988661749SPeter Brune     if (by) {
39088661749SPeter Brune       Ys = ys - s;
39188661749SPeter Brune     } else {
39288661749SPeter Brune       Ys = 0;
39388661749SPeter Brune     }
39488661749SPeter Brune     IYs = 0;
39588661749SPeter Brune   }
396d9c9ebe5SPeter Brune   if (ye+s <= N) {
397d9c9ebe5SPeter Brune     Ye = ye + s; IYe = ye + s;
39888661749SPeter Brune   } else {
39988661749SPeter Brune     if (by) {
40088661749SPeter Brune       Ye = ye + s;
40188661749SPeter Brune     } else {
40288661749SPeter Brune       Ye = N;
40388661749SPeter Brune     }
40488661749SPeter Brune     IYe = N;
40588661749SPeter Brune   }
40688661749SPeter Brune 
407bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) {
408d9c9ebe5SPeter Brune     IYs = ys - s;
409d9c9ebe5SPeter Brune     IYe = ye + s;
410d9c9ebe5SPeter Brune     Ys  = ys - s;
411d9c9ebe5SPeter Brune     Ye  = ye + s;
41288661749SPeter Brune   }
41388661749SPeter Brune 
41488661749SPeter Brune   /* stencil length in each direction */
415d9c9ebe5SPeter Brune   s_x = s;
416d9c9ebe5SPeter Brune   s_y = s;
41747c6ae99SBarry Smith 
41847c6ae99SBarry Smith   /* determine starting point of each processor */
41947c6ae99SBarry Smith   nn       = x*y;
420dcca6d9dSJed Brown   ierr     = PetscMalloc2(size+1,&bases,size,&ldims);CHKERRQ(ierr);
42147c6ae99SBarry Smith   ierr     = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
42247c6ae99SBarry Smith   bases[0] = 0;
42347c6ae99SBarry Smith   for (i=1; i<=size; i++) {
42447c6ae99SBarry Smith     bases[i] = ldims[i-1];
42547c6ae99SBarry Smith   }
42647c6ae99SBarry Smith   for (i=1; i<=size; i++) {
42747c6ae99SBarry Smith     bases[i] += bases[i-1];
42847c6ae99SBarry Smith   }
429ce00eea3SSatish Balay   base = bases[rank]*dof;
43047c6ae99SBarry Smith 
43147c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
432ce00eea3SSatish Balay   dd->Nlocal = x*y*dof;
433b1fb7eb7SBarry Smith   ierr       = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global);CHKERRQ(ierr);
434ce00eea3SSatish Balay   dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof;
435b1fb7eb7SBarry Smith   ierr       = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local);CHKERRQ(ierr);
43647c6ae99SBarry Smith 
43747c6ae99SBarry Smith   /* generate appropriate vector scatters */
43847c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
43902fe608eSBarry Smith   ierr  = PetscMalloc1((IXe-IXs)*(IYe-IYs),&idx);CHKERRQ(ierr);
440ce00eea3SSatish Balay   left  = xs - Xs; right = left + x;
441ce00eea3SSatish Balay   down  = ys - Ys; up = down + y;
44247c6ae99SBarry Smith   count = 0;
44347c6ae99SBarry Smith   for (i=down; i<up; i++) {
444ce00eea3SSatish Balay     for (j=left; j<right; j++) {
445ce00eea3SSatish Balay       idx[count++] = i*(Xe-Xs) + j;
44647c6ae99SBarry Smith     }
44747c6ae99SBarry Smith   }
44847c6ae99SBarry Smith 
449ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
450ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
451d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_BOX) {
452ce00eea3SSatish Balay     left  = IXs - Xs; right = left + (IXe-IXs);
453ce00eea3SSatish Balay     down  = IYs - Ys; up = down + (IYe-IYs);
454ce00eea3SSatish Balay     count = 0;
455ce00eea3SSatish Balay     for (i=down; i<up; i++) {
456ce00eea3SSatish Balay       for (j=left; j<right; j++) {
457ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
458ce00eea3SSatish Balay       }
459ce00eea3SSatish Balay     }
460ce00eea3SSatish Balay     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
461ce00eea3SSatish Balay 
46247c6ae99SBarry Smith   } else {
46347c6ae99SBarry Smith     /* must drop into cross shape region */
46447c6ae99SBarry Smith     /*       ---------|
46547c6ae99SBarry Smith             |  top    |
466ce00eea3SSatish Balay          |---         ---| up
46747c6ae99SBarry Smith          |   middle      |
46847c6ae99SBarry Smith          |               |
469ce00eea3SSatish Balay          ----         ---- down
47047c6ae99SBarry Smith             | bottom  |
47147c6ae99SBarry Smith             -----------
47247c6ae99SBarry Smith          Xs xs        xe Xe */
473ce00eea3SSatish Balay     left  = xs - Xs; right = left + x;
474ce00eea3SSatish Balay     down  = ys - Ys; up = down + y;
47547c6ae99SBarry Smith     count = 0;
476ce00eea3SSatish Balay     /* bottom */
477ce00eea3SSatish Balay     for (i=(IYs-Ys); i<down; i++) {
478ce00eea3SSatish Balay       for (j=left; j<right; j++) {
479ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
48047c6ae99SBarry Smith       }
48147c6ae99SBarry Smith     }
48247c6ae99SBarry Smith     /* middle */
48347c6ae99SBarry Smith     for (i=down; i<up; i++) {
484ce00eea3SSatish Balay       for (j=(IXs-Xs); j<(IXe-Xs); j++) {
485ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
48647c6ae99SBarry Smith       }
48747c6ae99SBarry Smith     }
48847c6ae99SBarry Smith     /* top */
489ce00eea3SSatish Balay     for (i=up; i<up+IYe-ye; i++) {
490ce00eea3SSatish Balay       for (j=left; j<right; j++) {
491ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
49247c6ae99SBarry Smith       }
49347c6ae99SBarry Smith     }
49447c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
49547c6ae99SBarry Smith   }
49647c6ae99SBarry Smith 
49747c6ae99SBarry Smith 
49847c6ae99SBarry Smith   /* determine who lies on each side of us stored in    n6 n7 n8
49947c6ae99SBarry Smith                                                         n3    n5
50047c6ae99SBarry Smith                                                         n0 n1 n2
50147c6ae99SBarry Smith   */
50247c6ae99SBarry Smith 
50347c6ae99SBarry Smith   /* Assume the Non-Periodic Case */
50447c6ae99SBarry Smith   n1 = rank - m;
50547c6ae99SBarry Smith   if (rank % m) {
50647c6ae99SBarry Smith     n0 = n1 - 1;
50747c6ae99SBarry Smith   } else {
50847c6ae99SBarry Smith     n0 = -1;
50947c6ae99SBarry Smith   }
51047c6ae99SBarry Smith   if ((rank+1) % m) {
51147c6ae99SBarry Smith     n2 = n1 + 1;
51247c6ae99SBarry Smith     n5 = rank + 1;
51347c6ae99SBarry Smith     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
51447c6ae99SBarry Smith   } else {
51547c6ae99SBarry Smith     n2 = -1; n5 = -1; n8 = -1;
51647c6ae99SBarry Smith   }
51747c6ae99SBarry Smith   if (rank % m) {
51847c6ae99SBarry Smith     n3 = rank - 1;
51947c6ae99SBarry Smith     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
52047c6ae99SBarry Smith   } else {
52147c6ae99SBarry Smith     n3 = -1; n6 = -1;
52247c6ae99SBarry Smith   }
52347c6ae99SBarry Smith   n7 = rank + m; if (n7 >= m*n) n7 = -1;
52447c6ae99SBarry Smith 
525bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC && by == DM_BOUNDARY_PERIODIC) {
52647c6ae99SBarry Smith     /* Modify for Periodic Cases */
52747c6ae99SBarry Smith     /* Handle all four corners */
52847c6ae99SBarry Smith     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
52947c6ae99SBarry Smith     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
53047c6ae99SBarry Smith     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
53147c6ae99SBarry Smith     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
53247c6ae99SBarry Smith 
53347c6ae99SBarry Smith     /* Handle Top and Bottom Sides */
53447c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
53547c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
53647c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
53747c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
53847c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
53947c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
54047c6ae99SBarry Smith 
54147c6ae99SBarry Smith     /* Handle Left and Right Sides */
54247c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
54347c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
54447c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
54547c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
54647c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
54747c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
548bff4a2f0SMatthew G. Knepley   } else if (by == DM_BOUNDARY_PERIODIC) {  /* Handle Top and Bottom Sides */
549ce00eea3SSatish Balay     if (n1 < 0) n1 = rank + m * (n-1);
550ce00eea3SSatish Balay     if (n7 < 0) n7 = rank - m * (n-1);
551ce00eea3SSatish Balay     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
552ce00eea3SSatish Balay     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
553ce00eea3SSatish Balay     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
554ce00eea3SSatish Balay     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
555bff4a2f0SMatthew G. Knepley   } else if (bx == DM_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */
556ce00eea3SSatish Balay     if (n3 < 0) n3 = rank + (m-1);
557ce00eea3SSatish Balay     if (n5 < 0) n5 = rank - (m-1);
558ce00eea3SSatish Balay     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
559ce00eea3SSatish Balay     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
560ce00eea3SSatish Balay     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
561ce00eea3SSatish Balay     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
56247c6ae99SBarry Smith   }
563ce00eea3SSatish Balay 
564785e854fSJed Brown   ierr = PetscMalloc1(9,&dd->neighbors);CHKERRQ(ierr);
5658865f1eaSKarl Rupp 
56647c6ae99SBarry Smith   dd->neighbors[0] = n0;
56747c6ae99SBarry Smith   dd->neighbors[1] = n1;
56847c6ae99SBarry Smith   dd->neighbors[2] = n2;
56947c6ae99SBarry Smith   dd->neighbors[3] = n3;
57047c6ae99SBarry Smith   dd->neighbors[4] = rank;
57147c6ae99SBarry Smith   dd->neighbors[5] = n5;
57247c6ae99SBarry Smith   dd->neighbors[6] = n6;
57347c6ae99SBarry Smith   dd->neighbors[7] = n7;
57447c6ae99SBarry Smith   dd->neighbors[8] = n8;
57547c6ae99SBarry Smith 
576d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
57747c6ae99SBarry Smith     /* save corner processor numbers */
57847c6ae99SBarry Smith     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
57947c6ae99SBarry Smith     n0  = n2 = n6 = n8 = -1;
58047c6ae99SBarry Smith   }
58147c6ae99SBarry Smith 
582785e854fSJed Brown   ierr = PetscMalloc1((Xe-Xs)*(Ye-Ys),&idx);CHKERRQ(ierr);
58347c6ae99SBarry Smith 
584ce00eea3SSatish Balay   nn = 0;
58547c6ae99SBarry Smith   xbase = bases[rank];
58647c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
58747c6ae99SBarry Smith     if (n0 >= 0) { /* left below */
588ce00eea3SSatish Balay       x_t = lx[n0 % m];
58947c6ae99SBarry Smith       y_t = ly[(n0/m)];
59047c6ae99SBarry Smith       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
5918865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
59247c6ae99SBarry Smith     }
593ac119b13SBarry Smith 
59447c6ae99SBarry Smith     if (n1 >= 0) { /* directly below */
59547c6ae99SBarry Smith       x_t = x;
59647c6ae99SBarry Smith       y_t = ly[(n1/m)];
59747c6ae99SBarry Smith       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
5988865f1eaSKarl Rupp       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
599bff4a2f0SMatthew G. Knepley     } else if (by == DM_BOUNDARY_MIRROR) {
6008865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
60147c6ae99SBarry Smith     }
602ac119b13SBarry Smith 
60347c6ae99SBarry Smith     if (n2 >= 0) { /* right below */
604ce00eea3SSatish Balay       x_t = lx[n2 % m];
60547c6ae99SBarry Smith       y_t = ly[(n2/m)];
60647c6ae99SBarry Smith       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
6078865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
60847c6ae99SBarry Smith     }
60947c6ae99SBarry Smith   }
61047c6ae99SBarry Smith 
61147c6ae99SBarry Smith   for (i=0; i<y; i++) {
61247c6ae99SBarry Smith     if (n3 >= 0) { /* directly left */
613ce00eea3SSatish Balay       x_t = lx[n3 % m];
61447c6ae99SBarry Smith       /* y_t = y; */
61547c6ae99SBarry Smith       s_t = bases[n3] + (i+1)*x_t - s_x;
6168865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
617bff4a2f0SMatthew G. Knepley     } else if (bx == DM_BOUNDARY_MIRROR) {
6188865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
61947c6ae99SBarry Smith     }
62047c6ae99SBarry Smith 
6218865f1eaSKarl Rupp     for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
62247c6ae99SBarry Smith 
62347c6ae99SBarry Smith     if (n5 >= 0) { /* directly right */
624ce00eea3SSatish Balay       x_t = lx[n5 % m];
62547c6ae99SBarry Smith       /* y_t = y; */
62647c6ae99SBarry Smith       s_t = bases[n5] + (i)*x_t;
6278865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
628bff4a2f0SMatthew G. Knepley     } else if (bx == DM_BOUNDARY_MIRROR) {
6298865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
63047c6ae99SBarry Smith     }
63147c6ae99SBarry Smith   }
63247c6ae99SBarry Smith 
63347c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
63447c6ae99SBarry Smith     if (n6 >= 0) { /* left above */
635ce00eea3SSatish Balay       x_t = lx[n6 % m];
63647c6ae99SBarry Smith       /* y_t = ly[(n6/m)]; */
63747c6ae99SBarry Smith       s_t = bases[n6] + (i)*x_t - s_x;
6388865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
63947c6ae99SBarry Smith     }
640ac119b13SBarry Smith 
64147c6ae99SBarry Smith     if (n7 >= 0) { /* directly above */
64247c6ae99SBarry Smith       x_t = x;
64347c6ae99SBarry Smith       /* y_t = ly[(n7/m)]; */
64447c6ae99SBarry Smith       s_t = bases[n7] + (i-1)*x_t;
6458865f1eaSKarl Rupp       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
646bff4a2f0SMatthew G. Knepley     } else if (by == DM_BOUNDARY_MIRROR) {
6478865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
64847c6ae99SBarry Smith     }
649ac119b13SBarry Smith 
65047c6ae99SBarry Smith     if (n8 >= 0) { /* right above */
651ce00eea3SSatish Balay       x_t = lx[n8 % m];
65247c6ae99SBarry Smith       /* y_t = ly[(n8/m)]; */
65347c6ae99SBarry Smith       s_t = bases[n8] + (i-1)*x_t;
6548865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
65547c6ae99SBarry Smith     }
65647c6ae99SBarry Smith   }
65747c6ae99SBarry Smith 
658b1fb7eb7SBarry Smith   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from);CHKERRQ(ierr);
65947c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
6603bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr);
661fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
662fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
66347c6ae99SBarry Smith 
664d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
665ce00eea3SSatish Balay     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
666ce00eea3SSatish Balay   }
667ce00eea3SSatish Balay 
66888661749SPeter Brune   if (((stencil_type == DMDA_STENCIL_STAR)  ||
669bff4a2f0SMatthew G. Knepley        (bx && bx != DM_BOUNDARY_PERIODIC) ||
670bff4a2f0SMatthew G. Knepley        (by && by != DM_BOUNDARY_PERIODIC))) {
67147c6ae99SBarry Smith     /*
67247c6ae99SBarry Smith         Recompute the local to global mappings, this time keeping the
673ce00eea3SSatish Balay       information about the cross corner processor numbers and any ghosted
674ce00eea3SSatish Balay       but not periodic indices.
67547c6ae99SBarry Smith     */
67647c6ae99SBarry Smith     nn    = 0;
67747c6ae99SBarry Smith     xbase = bases[rank];
67847c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
67947c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
680ce00eea3SSatish Balay         x_t = lx[n0 % m];
68147c6ae99SBarry Smith         y_t = ly[(n0/m)];
68247c6ae99SBarry Smith         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
6838865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
684ce00eea3SSatish Balay       } else if (xs-Xs > 0 && ys-Ys > 0) {
6858865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
68647c6ae99SBarry Smith       }
68747c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
68847c6ae99SBarry Smith         x_t = x;
68947c6ae99SBarry Smith         y_t = ly[(n1/m)];
69047c6ae99SBarry Smith         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
6918865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
692ce00eea3SSatish Balay       } else if (ys-Ys > 0) {
693bff4a2f0SMatthew G. Knepley         if (by == DM_BOUNDARY_MIRROR) {
6948865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
695624904c4SBarry Smith         } else {
6968865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
69747c6ae99SBarry Smith         }
698624904c4SBarry Smith       }
69947c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
700ce00eea3SSatish Balay         x_t = lx[n2 % m];
70147c6ae99SBarry Smith         y_t = ly[(n2/m)];
70247c6ae99SBarry Smith         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
7038865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
704ce00eea3SSatish Balay       } else if (Xe-xe> 0 && ys-Ys > 0) {
7058865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
70647c6ae99SBarry Smith       }
70747c6ae99SBarry Smith     }
70847c6ae99SBarry Smith 
70947c6ae99SBarry Smith     for (i=0; i<y; i++) {
71047c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
711ce00eea3SSatish Balay         x_t = lx[n3 % m];
71247c6ae99SBarry Smith         /* y_t = y; */
71347c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x;
7148865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
715ce00eea3SSatish Balay       } else if (xs-Xs > 0) {
716bff4a2f0SMatthew G. Knepley         if (bx == DM_BOUNDARY_MIRROR) {
7178865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
718624904c4SBarry Smith         } else {
7198865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
72047c6ae99SBarry Smith         }
721624904c4SBarry Smith       }
72247c6ae99SBarry Smith 
7238865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
72447c6ae99SBarry Smith 
72547c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
726ce00eea3SSatish Balay         x_t = lx[n5 % m];
72747c6ae99SBarry Smith         /* y_t = y; */
72847c6ae99SBarry Smith         s_t = bases[n5] + (i)*x_t;
7298865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
730ce00eea3SSatish Balay       } else if (Xe-xe > 0) {
731bff4a2f0SMatthew G. Knepley         if (bx == DM_BOUNDARY_MIRROR) {
7328865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
733624904c4SBarry Smith         } else {
7348865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
73547c6ae99SBarry Smith         }
73647c6ae99SBarry Smith       }
737624904c4SBarry Smith     }
73847c6ae99SBarry Smith 
73947c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
74047c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
741ce00eea3SSatish Balay         x_t = lx[n6 % m];
74247c6ae99SBarry Smith         /* y_t = ly[(n6/m)]; */
74347c6ae99SBarry Smith         s_t = bases[n6] + (i)*x_t - s_x;
7448865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
745ce00eea3SSatish Balay       } else if (xs-Xs > 0 && Ye-ye > 0) {
7468865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
74747c6ae99SBarry Smith       }
74847c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
74947c6ae99SBarry Smith         x_t = x;
75047c6ae99SBarry Smith         /* y_t = ly[(n7/m)]; */
75147c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t;
7528865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
753ce00eea3SSatish Balay       } else if (Ye-ye > 0) {
754bff4a2f0SMatthew G. Knepley         if (by == DM_BOUNDARY_MIRROR) {
7558865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
756624904c4SBarry Smith         } else {
7578865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
75847c6ae99SBarry Smith         }
759624904c4SBarry Smith       }
76047c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
761ce00eea3SSatish Balay         x_t = lx[n8 % m];
76247c6ae99SBarry Smith         /* y_t = ly[(n8/m)]; */
76347c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t;
7648865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
765ce00eea3SSatish Balay       } else if (Xe-xe > 0 && Ye-ye > 0) {
7668865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
76747c6ae99SBarry Smith       }
76847c6ae99SBarry Smith     }
76947c6ae99SBarry Smith   }
770ce00eea3SSatish Balay   /*
771ce00eea3SSatish Balay      Set the local to global ordering in the global vector, this allows use
772ce00eea3SSatish Balay      of VecSetValuesLocal().
773ce00eea3SSatish Balay   */
77445b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr);
7753bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr);
77647c6ae99SBarry Smith 
777ce00eea3SSatish Balay   ierr  = PetscFree2(bases,ldims);CHKERRQ(ierr);
77847c6ae99SBarry Smith   dd->m = m;  dd->n  = n;
779ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
780ce00eea3SSatish Balay   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
781ce00eea3SSatish Balay   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
78247c6ae99SBarry Smith 
783fcfd50ebSBarry Smith   ierr = VecDestroy(&local);CHKERRQ(ierr);
784fcfd50ebSBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
78547c6ae99SBarry Smith 
78647c6ae99SBarry Smith   dd->gtol      = gtol;
78747c6ae99SBarry Smith   dd->base      = base;
7889a42bb27SBarry Smith   da->ops->view = DMView_DA_2d;
7890298fd71SBarry Smith   dd->ltol      = NULL;
7900298fd71SBarry Smith   dd->ao        = NULL;
79147c6ae99SBarry Smith   PetscFunctionReturn(0);
79247c6ae99SBarry Smith }
79347c6ae99SBarry Smith 
79447c6ae99SBarry Smith #undef __FUNCT__
795aa219208SBarry Smith #define __FUNCT__ "DMDACreate2d"
79647c6ae99SBarry Smith /*@C
797aa219208SBarry Smith    DMDACreate2d -  Creates an object that will manage the communication of  two-dimensional
79847c6ae99SBarry Smith    regular array data that is distributed across some processors.
79947c6ae99SBarry Smith 
80047c6ae99SBarry Smith    Collective on MPI_Comm
80147c6ae99SBarry Smith 
80247c6ae99SBarry Smith    Input Parameters:
80347c6ae99SBarry Smith +  comm - MPI communicator
8041321219cSEthan Coon .  bx,by - type of ghost nodes the array have.
805bff4a2f0SMatthew G. Knepley          Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC.
806aa219208SBarry Smith .  stencil_type - stencil type.  Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
80747c6ae99SBarry 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
80847c6ae99SBarry Smith             from the command line with -da_grid_x <M> -da_grid_y <N>)
80947c6ae99SBarry Smith .  m,n - corresponding number of processors in each dimension
81047c6ae99SBarry Smith          (or PETSC_DECIDE to have calculated)
81147c6ae99SBarry Smith .  dof - number of degrees of freedom per node
81247c6ae99SBarry Smith .  s - stencil width
81347c6ae99SBarry Smith -  lx, ly - arrays containing the number of nodes in each cell along
8140298fd71SBarry Smith            the x and y coordinates, or NULL. If non-null, these
81547c6ae99SBarry Smith            must be of length as m and n, and the corresponding
81647c6ae99SBarry Smith            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
81747c6ae99SBarry Smith            must be M, and the sum of the ly[] entries must be N.
81847c6ae99SBarry Smith 
81947c6ae99SBarry Smith    Output Parameter:
82047c6ae99SBarry Smith .  da - the resulting distributed array object
82147c6ae99SBarry Smith 
82247c6ae99SBarry Smith    Options Database Key:
823706a11cbSBarry Smith +  -dm_view - Calls DMView() at the conclusion of DMDACreate2d()
82447c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
82547c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
82647c6ae99SBarry Smith .  -da_processors_x <nx> - number of processors in x direction
82747c6ae99SBarry Smith .  -da_processors_y <ny> - number of processors in y direction
828e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement ratio in x direction
829e0f5d30fSBarry Smith .  -da_refine_y <ry> - refinement ratio in y direction
830e0f5d30fSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating, if M or N < 0
831e0f5d30fSBarry Smith 
83247c6ae99SBarry Smith 
83347c6ae99SBarry Smith    Level: beginner
83447c6ae99SBarry Smith 
83547c6ae99SBarry Smith    Notes:
836aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
837aa219208SBarry Smith    standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
83847c6ae99SBarry Smith    the standard 9-pt stencil.
83947c6ae99SBarry Smith 
840aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
841564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
842564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
84347c6ae99SBarry Smith 
84447c6ae99SBarry Smith .keywords: distributed array, create, two-dimensional
84547c6ae99SBarry Smith 
846aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
84799f0b487SRichard Tran Mills           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(),
848d461ba97SBarry Smith           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
84947c6ae99SBarry Smith 
85047c6ae99SBarry Smith @*/
851fe16a2e9SBarry Smith 
852bff4a2f0SMatthew G. Knepley PetscErrorCode  DMDACreate2d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMDAStencilType stencil_type,
8539a42bb27SBarry Smith                           PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
85447c6ae99SBarry Smith {
85547c6ae99SBarry Smith   PetscErrorCode ierr;
85647c6ae99SBarry Smith 
85747c6ae99SBarry Smith   PetscFunctionBegin;
858aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
859c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*da, 2);CHKERRQ(ierr);
860aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, 1);CHKERRQ(ierr);
861aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr);
862bff4a2f0SMatthew G. Knepley   ierr = DMDASetBoundaryType(*da, bx, by, DM_BOUNDARY_NONE);CHKERRQ(ierr);
863aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
864aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
865aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
8660298fd71SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, NULL);CHKERRQ(ierr);
86747c6ae99SBarry Smith   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
8689a42bb27SBarry Smith   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
8699a42bb27SBarry Smith   ierr = DMSetUp(*da);CHKERRQ(ierr);
87047c6ae99SBarry Smith   PetscFunctionReturn(0);
87147c6ae99SBarry Smith }
872