xref: /petsc/src/dm/impls/da/da2.c (revision 8135c375e96a7bf70d67c02f02cd826154817056)
19a42bb27SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
39804daf3SBarry Smith #include <petscdraw.h>
447c6ae99SBarry Smith 
5e0877f53SBarry Smith static PetscErrorCode DMView_DA_2d(DM da,PetscViewer viewer)
647c6ae99SBarry Smith {
747c6ae99SBarry Smith   PetscErrorCode ierr;
847c6ae99SBarry Smith   PetscMPIInt    rank;
99a42bb27SBarry Smith   PetscBool      iascii,isdraw,isbinary;
1047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
119a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
129a42bb27SBarry Smith   PetscBool ismatlab;
139a42bb27SBarry Smith #endif
1447c6ae99SBarry Smith 
1547c6ae99SBarry Smith   PetscFunctionBegin;
16ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr);
1747c6ae99SBarry Smith 
18251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
19251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
20251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
219a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
22251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr);
239a42bb27SBarry Smith #endif
2447c6ae99SBarry Smith   if (iascii) {
2547c6ae99SBarry Smith     PetscViewerFormat format;
2647c6ae99SBarry Smith 
2747c6ae99SBarry Smith     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
28*8135c375SStefano Zampini     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL && format != PETSC_VIEWER_ASCII_GLVIS) {
29aa219208SBarry Smith       DMDALocalInfo info;
30aa219208SBarry Smith       ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
311575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
3247c6ae99SBarry 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);
3347c6ae99SBarry 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);
3447c6ae99SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
351575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
36*8135c375SStefano Zampini     } else if (format == PETSC_VIEWER_ASCII_GLVIS) {
37*8135c375SStefano Zampini       ierr = DMView_DA_GLVis(da,viewer);CHKERRQ(ierr);
383da9ae13SJed Brown     } else {
393da9ae13SJed Brown       ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr);
4047c6ae99SBarry Smith     }
4147c6ae99SBarry Smith   } else if (isdraw) {
4247c6ae99SBarry Smith     PetscDraw      draw;
4347c6ae99SBarry Smith     double         ymin = -1*dd->s-1,ymax = dd->N+dd->s;
4447c6ae99SBarry Smith     double         xmin = -1*dd->s-1,xmax = dd->M+dd->s;
4547c6ae99SBarry Smith     double         x,y;
468ea3bf28SBarry Smith     PetscInt       base;
478ea3bf28SBarry Smith     const PetscInt *idx;
4847c6ae99SBarry Smith     char           node[10];
4947c6ae99SBarry Smith     PetscBool      isnull;
5047c6ae99SBarry Smith 
5147c6ae99SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
525b399a63SLisandro Dalcin     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
535b399a63SLisandro Dalcin     if (isnull) PetscFunctionReturn(0);
5447c6ae99SBarry Smith 
555b399a63SLisandro Dalcin     ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
565b399a63SLisandro Dalcin     ierr = PetscDrawClear(draw);CHKERRQ(ierr);
575b399a63SLisandro Dalcin     ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
585b399a63SLisandro Dalcin 
595b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
6047c6ae99SBarry Smith     /* first processor draw all node lines */
6147c6ae99SBarry Smith     if (!rank) {
6247c6ae99SBarry Smith       ymin = 0.0; ymax = dd->N - 1;
6347c6ae99SBarry Smith       for (xmin=0; xmin<dd->M; xmin++) {
6447c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
6547c6ae99SBarry Smith       }
6647c6ae99SBarry Smith       xmin = 0.0; xmax = dd->M - 1;
6747c6ae99SBarry Smith       for (ymin=0; ymin<dd->N; ymin++) {
6847c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
6947c6ae99SBarry Smith       }
7047c6ae99SBarry Smith     }
715b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
725b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
7347c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
7447c6ae99SBarry Smith 
755b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
7647c6ae99SBarry Smith     /* draw my box */
775b399a63SLisandro Dalcin     xmin = dd->xs/dd->w; xmax =(dd->xe-1)/dd->w; ymin = dd->ys; ymax = dd->ye - 1;
7847c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
7947c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
8047c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
8147c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
8247c6ae99SBarry Smith     /* put in numbers */
8347c6ae99SBarry Smith     base = (dd->base)/dd->w;
8447c6ae99SBarry Smith     for (y=ymin; y<=ymax; y++) {
8547c6ae99SBarry Smith       for (x=xmin; x<=xmax; x++) {
865b399a63SLisandro Dalcin         ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)base++);CHKERRQ(ierr);
8747c6ae99SBarry Smith         ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
8847c6ae99SBarry Smith       }
8947c6ae99SBarry Smith     }
905b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
915b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
9247c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
9347c6ae99SBarry Smith 
945b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
955b399a63SLisandro Dalcin     /* overlay ghost numbers, useful for error checking */
9645b6f7e9SBarry Smith     ierr = ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr);
975b399a63SLisandro Dalcin     base = 0; xmin = dd->Xs; xmax = dd->Xe; ymin = dd->Ys; ymax = dd->Ye;
9847c6ae99SBarry Smith     for (y=ymin; y<ymax; y++) {
9947c6ae99SBarry Smith       for (x=xmin; x<xmax; x++) {
10047c6ae99SBarry Smith         if ((base % dd->w) == 0) {
1015b399a63SLisandro Dalcin           ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)(idx[base/dd->w]));CHKERRQ(ierr);
10247c6ae99SBarry Smith           ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
10347c6ae99SBarry Smith         }
10447c6ae99SBarry Smith         base++;
10547c6ae99SBarry Smith       }
10647c6ae99SBarry Smith     }
107302440fdSBarry Smith     ierr = ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr);
1085b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
1095b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
11047c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
111832b7cebSLisandro Dalcin     ierr = PetscDrawSave(draw);CHKERRQ(ierr);
1129a42bb27SBarry Smith   } else if (isbinary) {
1139a42bb27SBarry Smith     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
1149a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1159a42bb27SBarry Smith   } else if (ismatlab) {
1169a42bb27SBarry Smith     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
1179a42bb27SBarry Smith #endif
11811aeaf0aSBarry Smith   }
11947c6ae99SBarry Smith   PetscFunctionReturn(0);
12047c6ae99SBarry Smith }
12147c6ae99SBarry Smith 
12247c6ae99SBarry Smith /*
12347c6ae99SBarry Smith       M is number of grid points
12447c6ae99SBarry Smith       m is number of processors
12547c6ae99SBarry Smith 
12647c6ae99SBarry Smith */
1277087cfbeSBarry Smith PetscErrorCode  DMDASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm)
12847c6ae99SBarry Smith {
12947c6ae99SBarry Smith   PetscErrorCode ierr;
13047c6ae99SBarry Smith   PetscInt       m,n = 0,x = 0,y = 0;
13147c6ae99SBarry Smith   PetscMPIInt    size,csize,rank;
13247c6ae99SBarry Smith 
13347c6ae99SBarry Smith   PetscFunctionBegin;
13447c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
13547c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13647c6ae99SBarry Smith 
13747c6ae99SBarry Smith   csize = 4*size;
13847c6ae99SBarry Smith   do {
13947c6ae99SBarry 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);
14047c6ae99SBarry Smith     csize = csize/4;
14147c6ae99SBarry Smith 
142369cc0aeSBarry Smith     m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)csize)/((PetscReal)N)));
14347c6ae99SBarry Smith     if (!m) m = 1;
14447c6ae99SBarry Smith     while (m > 0) {
14547c6ae99SBarry Smith       n = csize/m;
14647c6ae99SBarry Smith       if (m*n == csize) break;
14747c6ae99SBarry Smith       m--;
14847c6ae99SBarry Smith     }
14947c6ae99SBarry Smith     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
15047c6ae99SBarry Smith 
15147c6ae99SBarry Smith     x = M/m + ((M % m) > ((csize-1) % m));
15247c6ae99SBarry Smith     y = (N + (csize-1)/m)/n;
15347c6ae99SBarry Smith   } while ((x < 4 || y < 4) && csize > 1);
15447c6ae99SBarry Smith   if (size != csize) {
15547c6ae99SBarry Smith     MPI_Group   entire_group,sub_group;
15647c6ae99SBarry Smith     PetscMPIInt i,*groupies;
15747c6ae99SBarry Smith 
15847c6ae99SBarry Smith     ierr = MPI_Comm_group(comm,&entire_group);CHKERRQ(ierr);
159785e854fSJed Brown     ierr = PetscMalloc1(csize,&groupies);CHKERRQ(ierr);
16047c6ae99SBarry Smith     for (i=0; i<csize; i++) {
16147c6ae99SBarry Smith       groupies[i] = (rank/csize)*csize + i;
16247c6ae99SBarry Smith     }
16347c6ae99SBarry Smith     ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);CHKERRQ(ierr);
16447c6ae99SBarry Smith     ierr = PetscFree(groupies);CHKERRQ(ierr);
16547c6ae99SBarry Smith     ierr = MPI_Comm_create(comm,sub_group,outcomm);CHKERRQ(ierr);
16647c6ae99SBarry Smith     ierr = MPI_Group_free(&entire_group);CHKERRQ(ierr);
16747c6ae99SBarry Smith     ierr = MPI_Group_free(&sub_group);CHKERRQ(ierr);
168aa219208SBarry Smith     ierr = PetscInfo1(0,"DMDASplitComm2d:Creating redundant coarse problems of size %d\n",csize);CHKERRQ(ierr);
16947c6ae99SBarry Smith   } else {
17047c6ae99SBarry Smith     *outcomm = comm;
17147c6ae99SBarry Smith   }
17247c6ae99SBarry Smith   PetscFunctionReturn(0);
17347c6ae99SBarry Smith }
17447c6ae99SBarry Smith 
17547c6ae99SBarry Smith #if defined(new)
17647c6ae99SBarry Smith /*
177aa219208SBarry Smith   DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
178aa219208SBarry Smith     function lives on a DMDA
17947c6ae99SBarry Smith 
18047c6ae99SBarry Smith         y ~= (F(u + ha) - F(u))/h,
18147c6ae99SBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
18247c6ae99SBarry Smith         u = current iterate
18347c6ae99SBarry Smith         h = difference interval
18447c6ae99SBarry Smith */
185aa219208SBarry Smith PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a)
18647c6ae99SBarry Smith {
18747c6ae99SBarry Smith   PetscScalar    h,*aa,*ww,v;
18847c6ae99SBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
18947c6ae99SBarry Smith   PetscErrorCode ierr;
19047c6ae99SBarry Smith   PetscInt       gI,nI;
19147c6ae99SBarry Smith   MatStencil     stencil;
192aa219208SBarry Smith   DMDALocalInfo  info;
19347c6ae99SBarry Smith 
19447c6ae99SBarry Smith   PetscFunctionBegin;
19547c6ae99SBarry Smith   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
19647c6ae99SBarry Smith   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
19747c6ae99SBarry Smith 
19847c6ae99SBarry Smith   ierr = VecGetArray(U,&ww);CHKERRQ(ierr);
19947c6ae99SBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
20047c6ae99SBarry Smith 
20147c6ae99SBarry Smith   nI = 0;
20247c6ae99SBarry Smith   h  = ww[gI];
20347c6ae99SBarry Smith   if (h == 0.0) h = 1.0;
20447c6ae99SBarry Smith   if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin;
20547c6ae99SBarry Smith   else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
20647c6ae99SBarry Smith   h *= epsilon;
20747c6ae99SBarry Smith 
20847c6ae99SBarry Smith   ww[gI] += h;
20947c6ae99SBarry Smith   ierr    = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
21047c6ae99SBarry Smith   aa[nI]  = (v - aa[nI])/h;
21147c6ae99SBarry Smith   ww[gI] -= h;
21247c6ae99SBarry Smith   nI++;
2138865f1eaSKarl Rupp 
21447c6ae99SBarry Smith   ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr);
21547c6ae99SBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
21647c6ae99SBarry Smith   PetscFunctionReturn(0);
21747c6ae99SBarry Smith }
21847c6ae99SBarry Smith #endif
21947c6ae99SBarry Smith 
2207087cfbeSBarry Smith PetscErrorCode  DMSetUp_DA_2D(DM da)
22147c6ae99SBarry Smith {
22247c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
22347c6ae99SBarry Smith   const PetscInt   M            = dd->M;
22447c6ae99SBarry Smith   const PetscInt   N            = dd->N;
22547c6ae99SBarry Smith   PetscInt         m            = dd->m;
22647c6ae99SBarry Smith   PetscInt         n            = dd->n;
22747c6ae99SBarry Smith   const PetscInt   dof          = dd->w;
22847c6ae99SBarry Smith   const PetscInt   s            = dd->s;
229bff4a2f0SMatthew G. Knepley   DMBoundaryType   bx           = dd->bx;
230bff4a2f0SMatthew G. Knepley   DMBoundaryType   by           = dd->by;
23119fd82e9SBarry Smith   DMDAStencilType  stencil_type = dd->stencil_type;
23247c6ae99SBarry Smith   PetscInt         *lx          = dd->lx;
23347c6ae99SBarry Smith   PetscInt         *ly          = dd->ly;
23447c6ae99SBarry Smith   MPI_Comm         comm;
23547c6ae99SBarry Smith   PetscMPIInt      rank,size;
236bd1fc5aeSBarry Smith   PetscInt         xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,IXs,IXe,IYs,IYe;
2378ea3bf28SBarry Smith   PetscInt         up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn;
238db87c5ecSEthan Coon   PetscInt         xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
23947c6ae99SBarry Smith   PetscInt         s_x,s_y; /* s proportionalized to w */
24047c6ae99SBarry Smith   PetscInt         sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
24147c6ae99SBarry Smith   Vec              local,global;
242bd1fc5aeSBarry Smith   VecScatter       gtol;
24345b6f7e9SBarry Smith   IS               to,from;
24447c6ae99SBarry Smith   PetscErrorCode   ierr;
24547c6ae99SBarry Smith 
24647c6ae99SBarry Smith   PetscFunctionBegin;
247bff4a2f0SMatthew 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");
24847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
2493855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES)
250bafee8b4SSatish Balay   if (((PetscInt64) M)*((PetscInt64) N)*((PetscInt64) dof) > (PetscInt64) 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);
2513855c12bSBarry Smith #endif
25247c6ae99SBarry Smith 
25347c6ae99SBarry Smith   if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
25447c6ae99SBarry Smith   if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
25547c6ae99SBarry Smith 
25647c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
25747c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
25847c6ae99SBarry Smith 
2597d310018SBarry Smith   dd->p = 1;
26047c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
26147c6ae99SBarry Smith     if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
26247c6ae99SBarry Smith     else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
26347c6ae99SBarry Smith   }
26447c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
26547c6ae99SBarry Smith     if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
26647c6ae99SBarry Smith     else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
26747c6ae99SBarry Smith   }
26847c6ae99SBarry Smith 
26947c6ae99SBarry Smith   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
27047c6ae99SBarry Smith     if (n != PETSC_DECIDE) {
27147c6ae99SBarry Smith       m = size/n;
27247c6ae99SBarry Smith     } else if (m != PETSC_DECIDE) {
27347c6ae99SBarry Smith       n = size/m;
27447c6ae99SBarry Smith     } else {
27547c6ae99SBarry Smith       /* try for squarish distribution */
276369cc0aeSBarry Smith       m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N)));
27747c6ae99SBarry Smith       if (!m) m = 1;
27847c6ae99SBarry Smith       while (m > 0) {
27947c6ae99SBarry Smith         n = size/m;
28047c6ae99SBarry Smith         if (m*n == size) break;
28147c6ae99SBarry Smith         m--;
28247c6ae99SBarry Smith       }
28347c6ae99SBarry Smith       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
28447c6ae99SBarry Smith     }
28547c6ae99SBarry 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 ");
28647c6ae99SBarry Smith   } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
28747c6ae99SBarry Smith 
28847c6ae99SBarry Smith   if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
28947c6ae99SBarry Smith   if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
29047c6ae99SBarry Smith 
29147c6ae99SBarry Smith   /*
29247c6ae99SBarry Smith      Determine locally owned region
29347c6ae99SBarry Smith      xs is the first local node number, x is the number of local nodes
29447c6ae99SBarry Smith   */
29547c6ae99SBarry Smith   if (!lx) {
296785e854fSJed Brown     ierr = PetscMalloc1(m, &dd->lx);CHKERRQ(ierr);
29747c6ae99SBarry Smith     lx   = dd->lx;
29847c6ae99SBarry Smith     for (i=0; i<m; i++) {
29947c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > i);
30047c6ae99SBarry Smith     }
30147c6ae99SBarry Smith   }
30247c6ae99SBarry Smith   x  = lx[rank % m];
30347c6ae99SBarry Smith   xs = 0;
30447c6ae99SBarry Smith   for (i=0; i<(rank % m); i++) {
30547c6ae99SBarry Smith     xs += lx[i];
30647c6ae99SBarry Smith   }
30747c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
30847c6ae99SBarry Smith   left = xs;
30947c6ae99SBarry Smith   for (i=(rank % m); i<m; i++) {
31047c6ae99SBarry Smith     left += lx[i];
31147c6ae99SBarry Smith   }
31247c6ae99SBarry 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);
31347c6ae99SBarry Smith #endif
31447c6ae99SBarry Smith 
31547c6ae99SBarry Smith   /*
31647c6ae99SBarry Smith      Determine locally owned region
31747c6ae99SBarry Smith      ys is the first local node number, y is the number of local nodes
31847c6ae99SBarry Smith   */
31947c6ae99SBarry Smith   if (!ly) {
320785e854fSJed Brown     ierr = PetscMalloc1(n, &dd->ly);CHKERRQ(ierr);
32147c6ae99SBarry Smith     ly   = dd->ly;
32247c6ae99SBarry Smith     for (i=0; i<n; i++) {
32347c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > i);
32447c6ae99SBarry Smith     }
32547c6ae99SBarry Smith   }
32647c6ae99SBarry Smith   y  = ly[rank/m];
32747c6ae99SBarry Smith   ys = 0;
32847c6ae99SBarry Smith   for (i=0; i<(rank/m); i++) {
32947c6ae99SBarry Smith     ys += ly[i];
33047c6ae99SBarry Smith   }
33147c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
33247c6ae99SBarry Smith   left = ys;
33347c6ae99SBarry Smith   for (i=(rank/m); i<n; i++) {
33447c6ae99SBarry Smith     left += ly[i];
33547c6ae99SBarry Smith   }
33647c6ae99SBarry 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);
33747c6ae99SBarry Smith #endif
33847c6ae99SBarry Smith 
339bcea557cSEthan Coon   /*
340bcea557cSEthan Coon    check if the scatter requires more than one process neighbor or wraps around
341bcea557cSEthan Coon    the domain more than once
342bcea557cSEthan Coon   */
343bff4a2f0SMatthew 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);
344bff4a2f0SMatthew 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);
34547c6ae99SBarry Smith   xe = xs + x;
34647c6ae99SBarry Smith   ye = ys + y;
34747c6ae99SBarry Smith 
348ce00eea3SSatish Balay   /* determine ghost region (Xs) and region scattered into (IXs)  */
349d9c9ebe5SPeter Brune   if (xs-s > 0) {
350d9c9ebe5SPeter Brune     Xs = xs - s; IXs = xs - s;
35188661749SPeter Brune   } else {
35288661749SPeter Brune     if (bx) {
35388661749SPeter Brune       Xs = xs - s;
35488661749SPeter Brune     } else {
35588661749SPeter Brune       Xs = 0;
35688661749SPeter Brune     }
35788661749SPeter Brune     IXs = 0;
35888661749SPeter Brune   }
359d9c9ebe5SPeter Brune   if (xe+s <= M) {
360d9c9ebe5SPeter Brune     Xe = xe + s; IXe = xe + s;
36188661749SPeter Brune   } else {
36288661749SPeter Brune     if (bx) {
363d9c9ebe5SPeter Brune       Xs = xs - s; Xe = xe + s;
36488661749SPeter Brune     } else {
36588661749SPeter Brune       Xe = M;
36688661749SPeter Brune     }
36788661749SPeter Brune     IXe = M;
36888661749SPeter Brune   }
36947c6ae99SBarry Smith 
370bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) {
371d9c9ebe5SPeter Brune     IXs = xs - s;
372d9c9ebe5SPeter Brune     IXe = xe + s;
373d9c9ebe5SPeter Brune     Xs  = xs - s;
374d9c9ebe5SPeter Brune     Xe  = xe + s;
37588661749SPeter Brune   }
37647c6ae99SBarry Smith 
377d9c9ebe5SPeter Brune   if (ys-s > 0) {
378d9c9ebe5SPeter Brune     Ys = ys - s; IYs = ys - s;
37988661749SPeter Brune   } else {
38088661749SPeter Brune     if (by) {
38188661749SPeter Brune       Ys = ys - s;
38288661749SPeter Brune     } else {
38388661749SPeter Brune       Ys = 0;
38488661749SPeter Brune     }
38588661749SPeter Brune     IYs = 0;
38688661749SPeter Brune   }
387d9c9ebe5SPeter Brune   if (ye+s <= N) {
388d9c9ebe5SPeter Brune     Ye = ye + s; IYe = ye + s;
38988661749SPeter Brune   } else {
39088661749SPeter Brune     if (by) {
39188661749SPeter Brune       Ye = ye + s;
39288661749SPeter Brune     } else {
39388661749SPeter Brune       Ye = N;
39488661749SPeter Brune     }
39588661749SPeter Brune     IYe = N;
39688661749SPeter Brune   }
39788661749SPeter Brune 
398bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) {
399d9c9ebe5SPeter Brune     IYs = ys - s;
400d9c9ebe5SPeter Brune     IYe = ye + s;
401d9c9ebe5SPeter Brune     Ys  = ys - s;
402d9c9ebe5SPeter Brune     Ye  = ye + s;
40388661749SPeter Brune   }
40488661749SPeter Brune 
40588661749SPeter Brune   /* stencil length in each direction */
406d9c9ebe5SPeter Brune   s_x = s;
407d9c9ebe5SPeter Brune   s_y = s;
40847c6ae99SBarry Smith 
40947c6ae99SBarry Smith   /* determine starting point of each processor */
41047c6ae99SBarry Smith   nn       = x*y;
411dcca6d9dSJed Brown   ierr     = PetscMalloc2(size+1,&bases,size,&ldims);CHKERRQ(ierr);
41247c6ae99SBarry Smith   ierr     = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
41347c6ae99SBarry Smith   bases[0] = 0;
41447c6ae99SBarry Smith   for (i=1; i<=size; i++) {
41547c6ae99SBarry Smith     bases[i] = ldims[i-1];
41647c6ae99SBarry Smith   }
41747c6ae99SBarry Smith   for (i=1; i<=size; i++) {
41847c6ae99SBarry Smith     bases[i] += bases[i-1];
41947c6ae99SBarry Smith   }
420ce00eea3SSatish Balay   base = bases[rank]*dof;
42147c6ae99SBarry Smith 
42247c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
423ce00eea3SSatish Balay   dd->Nlocal = x*y*dof;
424b1fb7eb7SBarry Smith   ierr       = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global);CHKERRQ(ierr);
425ce00eea3SSatish Balay   dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof;
426b1fb7eb7SBarry Smith   ierr       = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local);CHKERRQ(ierr);
42747c6ae99SBarry Smith 
42847c6ae99SBarry Smith   /* generate appropriate vector scatters */
42947c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
43002fe608eSBarry Smith   ierr  = PetscMalloc1((IXe-IXs)*(IYe-IYs),&idx);CHKERRQ(ierr);
431ce00eea3SSatish Balay   left  = xs - Xs; right = left + x;
432ce00eea3SSatish Balay   down  = ys - Ys; up = down + y;
43347c6ae99SBarry Smith   count = 0;
43447c6ae99SBarry Smith   for (i=down; i<up; i++) {
435ce00eea3SSatish Balay     for (j=left; j<right; j++) {
436ce00eea3SSatish Balay       idx[count++] = i*(Xe-Xs) + j;
43747c6ae99SBarry Smith     }
43847c6ae99SBarry Smith   }
43947c6ae99SBarry Smith 
440ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
441ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
442d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_BOX) {
443ce00eea3SSatish Balay     left  = IXs - Xs; right = left + (IXe-IXs);
444ce00eea3SSatish Balay     down  = IYs - Ys; up = down + (IYe-IYs);
445ce00eea3SSatish Balay     count = 0;
446ce00eea3SSatish Balay     for (i=down; i<up; i++) {
447ce00eea3SSatish Balay       for (j=left; j<right; j++) {
448ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
449ce00eea3SSatish Balay       }
450ce00eea3SSatish Balay     }
451ce00eea3SSatish Balay     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
452ce00eea3SSatish Balay 
45347c6ae99SBarry Smith   } else {
45447c6ae99SBarry Smith     /* must drop into cross shape region */
45547c6ae99SBarry Smith     /*       ---------|
45647c6ae99SBarry Smith             |  top    |
457ce00eea3SSatish Balay          |---         ---| up
45847c6ae99SBarry Smith          |   middle      |
45947c6ae99SBarry Smith          |               |
460ce00eea3SSatish Balay          ----         ---- down
46147c6ae99SBarry Smith             | bottom  |
46247c6ae99SBarry Smith             -----------
46347c6ae99SBarry Smith          Xs xs        xe Xe */
464ce00eea3SSatish Balay     left  = xs - Xs; right = left + x;
465ce00eea3SSatish Balay     down  = ys - Ys; up = down + y;
46647c6ae99SBarry Smith     count = 0;
467ce00eea3SSatish Balay     /* bottom */
468ce00eea3SSatish Balay     for (i=(IYs-Ys); i<down; i++) {
469ce00eea3SSatish Balay       for (j=left; j<right; j++) {
470ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
47147c6ae99SBarry Smith       }
47247c6ae99SBarry Smith     }
47347c6ae99SBarry Smith     /* middle */
47447c6ae99SBarry Smith     for (i=down; i<up; i++) {
475ce00eea3SSatish Balay       for (j=(IXs-Xs); j<(IXe-Xs); j++) {
476ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
47747c6ae99SBarry Smith       }
47847c6ae99SBarry Smith     }
47947c6ae99SBarry Smith     /* top */
480ce00eea3SSatish Balay     for (i=up; i<up+IYe-ye; i++) {
481ce00eea3SSatish Balay       for (j=left; j<right; j++) {
482ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
48347c6ae99SBarry Smith       }
48447c6ae99SBarry Smith     }
48547c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
48647c6ae99SBarry Smith   }
48747c6ae99SBarry Smith 
48847c6ae99SBarry Smith 
48947c6ae99SBarry Smith   /* determine who lies on each side of us stored in    n6 n7 n8
49047c6ae99SBarry Smith                                                         n3    n5
49147c6ae99SBarry Smith                                                         n0 n1 n2
49247c6ae99SBarry Smith   */
49347c6ae99SBarry Smith 
49447c6ae99SBarry Smith   /* Assume the Non-Periodic Case */
49547c6ae99SBarry Smith   n1 = rank - m;
49647c6ae99SBarry Smith   if (rank % m) {
49747c6ae99SBarry Smith     n0 = n1 - 1;
49847c6ae99SBarry Smith   } else {
49947c6ae99SBarry Smith     n0 = -1;
50047c6ae99SBarry Smith   }
50147c6ae99SBarry Smith   if ((rank+1) % m) {
50247c6ae99SBarry Smith     n2 = n1 + 1;
50347c6ae99SBarry Smith     n5 = rank + 1;
50447c6ae99SBarry Smith     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
50547c6ae99SBarry Smith   } else {
50647c6ae99SBarry Smith     n2 = -1; n5 = -1; n8 = -1;
50747c6ae99SBarry Smith   }
50847c6ae99SBarry Smith   if (rank % m) {
50947c6ae99SBarry Smith     n3 = rank - 1;
51047c6ae99SBarry Smith     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
51147c6ae99SBarry Smith   } else {
51247c6ae99SBarry Smith     n3 = -1; n6 = -1;
51347c6ae99SBarry Smith   }
51447c6ae99SBarry Smith   n7 = rank + m; if (n7 >= m*n) n7 = -1;
51547c6ae99SBarry Smith 
516bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC && by == DM_BOUNDARY_PERIODIC) {
51747c6ae99SBarry Smith     /* Modify for Periodic Cases */
51847c6ae99SBarry Smith     /* Handle all four corners */
51947c6ae99SBarry Smith     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
52047c6ae99SBarry Smith     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
52147c6ae99SBarry Smith     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
52247c6ae99SBarry Smith     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
52347c6ae99SBarry Smith 
52447c6ae99SBarry Smith     /* Handle Top and Bottom Sides */
52547c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
52647c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
52747c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
52847c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
52947c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
53047c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
53147c6ae99SBarry Smith 
53247c6ae99SBarry Smith     /* Handle Left and Right Sides */
53347c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
53447c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
53547c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
53647c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
53747c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
53847c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
539bff4a2f0SMatthew G. Knepley   } else if (by == DM_BOUNDARY_PERIODIC) {  /* Handle Top and Bottom Sides */
540ce00eea3SSatish Balay     if (n1 < 0) n1 = rank + m * (n-1);
541ce00eea3SSatish Balay     if (n7 < 0) n7 = rank - m * (n-1);
542ce00eea3SSatish Balay     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
543ce00eea3SSatish Balay     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
544ce00eea3SSatish Balay     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
545ce00eea3SSatish Balay     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
546bff4a2f0SMatthew G. Knepley   } else if (bx == DM_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */
547ce00eea3SSatish Balay     if (n3 < 0) n3 = rank + (m-1);
548ce00eea3SSatish Balay     if (n5 < 0) n5 = rank - (m-1);
549ce00eea3SSatish Balay     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
550ce00eea3SSatish Balay     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
551ce00eea3SSatish Balay     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
552ce00eea3SSatish Balay     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
55347c6ae99SBarry Smith   }
554ce00eea3SSatish Balay 
555785e854fSJed Brown   ierr = PetscMalloc1(9,&dd->neighbors);CHKERRQ(ierr);
5568865f1eaSKarl Rupp 
55747c6ae99SBarry Smith   dd->neighbors[0] = n0;
55847c6ae99SBarry Smith   dd->neighbors[1] = n1;
55947c6ae99SBarry Smith   dd->neighbors[2] = n2;
56047c6ae99SBarry Smith   dd->neighbors[3] = n3;
56147c6ae99SBarry Smith   dd->neighbors[4] = rank;
56247c6ae99SBarry Smith   dd->neighbors[5] = n5;
56347c6ae99SBarry Smith   dd->neighbors[6] = n6;
56447c6ae99SBarry Smith   dd->neighbors[7] = n7;
56547c6ae99SBarry Smith   dd->neighbors[8] = n8;
56647c6ae99SBarry Smith 
567d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
56847c6ae99SBarry Smith     /* save corner processor numbers */
56947c6ae99SBarry Smith     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
57047c6ae99SBarry Smith     n0  = n2 = n6 = n8 = -1;
57147c6ae99SBarry Smith   }
57247c6ae99SBarry Smith 
573785e854fSJed Brown   ierr = PetscMalloc1((Xe-Xs)*(Ye-Ys),&idx);CHKERRQ(ierr);
57447c6ae99SBarry Smith 
575ce00eea3SSatish Balay   nn = 0;
57647c6ae99SBarry Smith   xbase = bases[rank];
57747c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
57847c6ae99SBarry Smith     if (n0 >= 0) { /* left below */
579ce00eea3SSatish Balay       x_t = lx[n0 % m];
58047c6ae99SBarry Smith       y_t = ly[(n0/m)];
58147c6ae99SBarry Smith       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
5828865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
58347c6ae99SBarry Smith     }
584ac119b13SBarry Smith 
58547c6ae99SBarry Smith     if (n1 >= 0) { /* directly below */
58647c6ae99SBarry Smith       x_t = x;
58747c6ae99SBarry Smith       y_t = ly[(n1/m)];
58847c6ae99SBarry Smith       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
5898865f1eaSKarl Rupp       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
590bff4a2f0SMatthew G. Knepley     } else if (by == DM_BOUNDARY_MIRROR) {
5918865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
59247c6ae99SBarry Smith     }
593ac119b13SBarry Smith 
59447c6ae99SBarry Smith     if (n2 >= 0) { /* right below */
595ce00eea3SSatish Balay       x_t = lx[n2 % m];
59647c6ae99SBarry Smith       y_t = ly[(n2/m)];
59747c6ae99SBarry Smith       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
5988865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
59947c6ae99SBarry Smith     }
60047c6ae99SBarry Smith   }
60147c6ae99SBarry Smith 
60247c6ae99SBarry Smith   for (i=0; i<y; i++) {
60347c6ae99SBarry Smith     if (n3 >= 0) { /* directly left */
604ce00eea3SSatish Balay       x_t = lx[n3 % m];
60547c6ae99SBarry Smith       /* y_t = y; */
60647c6ae99SBarry Smith       s_t = bases[n3] + (i+1)*x_t - s_x;
6078865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
608bff4a2f0SMatthew G. Knepley     } else if (bx == DM_BOUNDARY_MIRROR) {
6098865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
61047c6ae99SBarry Smith     }
61147c6ae99SBarry Smith 
6128865f1eaSKarl Rupp     for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
61347c6ae99SBarry Smith 
61447c6ae99SBarry Smith     if (n5 >= 0) { /* directly right */
615ce00eea3SSatish Balay       x_t = lx[n5 % m];
61647c6ae99SBarry Smith       /* y_t = y; */
61747c6ae99SBarry Smith       s_t = bases[n5] + (i)*x_t;
6188865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
619bff4a2f0SMatthew G. Knepley     } else if (bx == DM_BOUNDARY_MIRROR) {
6208865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
62147c6ae99SBarry Smith     }
62247c6ae99SBarry Smith   }
62347c6ae99SBarry Smith 
62447c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
62547c6ae99SBarry Smith     if (n6 >= 0) { /* left above */
626ce00eea3SSatish Balay       x_t = lx[n6 % m];
62747c6ae99SBarry Smith       /* y_t = ly[(n6/m)]; */
62847c6ae99SBarry Smith       s_t = bases[n6] + (i)*x_t - s_x;
6298865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
63047c6ae99SBarry Smith     }
631ac119b13SBarry Smith 
63247c6ae99SBarry Smith     if (n7 >= 0) { /* directly above */
63347c6ae99SBarry Smith       x_t = x;
63447c6ae99SBarry Smith       /* y_t = ly[(n7/m)]; */
63547c6ae99SBarry Smith       s_t = bases[n7] + (i-1)*x_t;
6368865f1eaSKarl Rupp       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
637bff4a2f0SMatthew G. Knepley     } else if (by == DM_BOUNDARY_MIRROR) {
6388865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
63947c6ae99SBarry Smith     }
640ac119b13SBarry Smith 
64147c6ae99SBarry Smith     if (n8 >= 0) { /* right above */
642ce00eea3SSatish Balay       x_t = lx[n8 % m];
64347c6ae99SBarry Smith       /* y_t = ly[(n8/m)]; */
64447c6ae99SBarry Smith       s_t = bases[n8] + (i-1)*x_t;
6458865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
64647c6ae99SBarry Smith     }
64747c6ae99SBarry Smith   }
64847c6ae99SBarry Smith 
649b1fb7eb7SBarry Smith   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from);CHKERRQ(ierr);
65047c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
6513bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr);
652fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
653fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
65447c6ae99SBarry Smith 
655d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
656ce00eea3SSatish Balay     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
657ce00eea3SSatish Balay   }
658ce00eea3SSatish Balay 
659288e7d53SBarry Smith   if (((stencil_type == DMDA_STENCIL_STAR)  || (bx && bx != DM_BOUNDARY_PERIODIC) || (by && by != DM_BOUNDARY_PERIODIC))) {
66047c6ae99SBarry Smith     /*
66147c6ae99SBarry Smith         Recompute the local to global mappings, this time keeping the
662ce00eea3SSatish Balay       information about the cross corner processor numbers and any ghosted
663ce00eea3SSatish Balay       but not periodic indices.
66447c6ae99SBarry Smith     */
66547c6ae99SBarry Smith     nn    = 0;
66647c6ae99SBarry Smith     xbase = bases[rank];
66747c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
66847c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
669ce00eea3SSatish Balay         x_t = lx[n0 % m];
67047c6ae99SBarry Smith         y_t = ly[(n0/m)];
67147c6ae99SBarry Smith         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
6728865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
673ce00eea3SSatish Balay       } else if (xs-Xs > 0 && ys-Ys > 0) {
6748865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
67547c6ae99SBarry Smith       }
67647c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
67747c6ae99SBarry Smith         x_t = x;
67847c6ae99SBarry Smith         y_t = ly[(n1/m)];
67947c6ae99SBarry Smith         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
6808865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
681ce00eea3SSatish Balay       } else if (ys-Ys > 0) {
682bff4a2f0SMatthew G. Knepley         if (by == DM_BOUNDARY_MIRROR) {
6838865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
684624904c4SBarry Smith         } else {
6858865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
68647c6ae99SBarry Smith         }
687624904c4SBarry Smith       }
68847c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
689ce00eea3SSatish Balay         x_t = lx[n2 % m];
69047c6ae99SBarry Smith         y_t = ly[(n2/m)];
69147c6ae99SBarry Smith         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
6928865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
693ce00eea3SSatish Balay       } else if (Xe-xe> 0 && ys-Ys > 0) {
6948865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
69547c6ae99SBarry Smith       }
69647c6ae99SBarry Smith     }
69747c6ae99SBarry Smith 
69847c6ae99SBarry Smith     for (i=0; i<y; i++) {
69947c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
700ce00eea3SSatish Balay         x_t = lx[n3 % m];
70147c6ae99SBarry Smith         /* y_t = y; */
70247c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x;
7038865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
704ce00eea3SSatish Balay       } else if (xs-Xs > 0) {
705bff4a2f0SMatthew G. Knepley         if (bx == DM_BOUNDARY_MIRROR) {
7068865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
707624904c4SBarry Smith         } else {
7088865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
70947c6ae99SBarry Smith         }
710624904c4SBarry Smith       }
71147c6ae99SBarry Smith 
7128865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
71347c6ae99SBarry Smith 
71447c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
715ce00eea3SSatish Balay         x_t = lx[n5 % m];
71647c6ae99SBarry Smith         /* y_t = y; */
71747c6ae99SBarry Smith         s_t = bases[n5] + (i)*x_t;
7188865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
719ce00eea3SSatish Balay       } else if (Xe-xe > 0) {
720bff4a2f0SMatthew G. Knepley         if (bx == DM_BOUNDARY_MIRROR) {
7218865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
722624904c4SBarry Smith         } else {
7238865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
72447c6ae99SBarry Smith         }
72547c6ae99SBarry Smith       }
726624904c4SBarry Smith     }
72747c6ae99SBarry Smith 
72847c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
72947c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
730ce00eea3SSatish Balay         x_t = lx[n6 % m];
73147c6ae99SBarry Smith         /* y_t = ly[(n6/m)]; */
73247c6ae99SBarry Smith         s_t = bases[n6] + (i)*x_t - s_x;
7338865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
734ce00eea3SSatish Balay       } else if (xs-Xs > 0 && Ye-ye > 0) {
7358865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
73647c6ae99SBarry Smith       }
73747c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
73847c6ae99SBarry Smith         x_t = x;
73947c6ae99SBarry Smith         /* y_t = ly[(n7/m)]; */
74047c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t;
7418865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
742ce00eea3SSatish Balay       } else if (Ye-ye > 0) {
743bff4a2f0SMatthew G. Knepley         if (by == DM_BOUNDARY_MIRROR) {
7448865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
745624904c4SBarry Smith         } else {
7468865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
74747c6ae99SBarry Smith         }
748624904c4SBarry Smith       }
74947c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
750ce00eea3SSatish Balay         x_t = lx[n8 % m];
75147c6ae99SBarry Smith         /* y_t = ly[(n8/m)]; */
75247c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t;
7538865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
754ce00eea3SSatish Balay       } else if (Xe-xe > 0 && Ye-ye > 0) {
7558865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
75647c6ae99SBarry Smith       }
75747c6ae99SBarry Smith     }
75847c6ae99SBarry Smith   }
759ce00eea3SSatish Balay   /*
760ce00eea3SSatish Balay      Set the local to global ordering in the global vector, this allows use
761ce00eea3SSatish Balay      of VecSetValuesLocal().
762ce00eea3SSatish Balay   */
76345b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr);
7643bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr);
76547c6ae99SBarry Smith 
766ce00eea3SSatish Balay   ierr  = PetscFree2(bases,ldims);CHKERRQ(ierr);
76747c6ae99SBarry Smith   dd->m = m;  dd->n  = n;
768ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
769ce00eea3SSatish Balay   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
770ce00eea3SSatish Balay   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
77147c6ae99SBarry Smith 
772fcfd50ebSBarry Smith   ierr = VecDestroy(&local);CHKERRQ(ierr);
773fcfd50ebSBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
77447c6ae99SBarry Smith 
77547c6ae99SBarry Smith   dd->gtol      = gtol;
77647c6ae99SBarry Smith   dd->base      = base;
7779a42bb27SBarry Smith   da->ops->view = DMView_DA_2d;
7780298fd71SBarry Smith   dd->ltol      = NULL;
7790298fd71SBarry Smith   dd->ao        = NULL;
78047c6ae99SBarry Smith   PetscFunctionReturn(0);
78147c6ae99SBarry Smith }
78247c6ae99SBarry Smith 
78347c6ae99SBarry Smith /*@C
784aa219208SBarry Smith    DMDACreate2d -  Creates an object that will manage the communication of  two-dimensional
78547c6ae99SBarry Smith    regular array data that is distributed across some processors.
78647c6ae99SBarry Smith 
78747c6ae99SBarry Smith    Collective on MPI_Comm
78847c6ae99SBarry Smith 
78947c6ae99SBarry Smith    Input Parameters:
79047c6ae99SBarry Smith +  comm - MPI communicator
7911321219cSEthan Coon .  bx,by - type of ghost nodes the array have.
792bff4a2f0SMatthew G. Knepley          Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC.
793aa219208SBarry Smith .  stencil_type - stencil type.  Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
794897f7067SBarry Smith .  M,N - global dimension in each direction of the array
79547c6ae99SBarry Smith .  m,n - corresponding number of processors in each dimension
79647c6ae99SBarry Smith          (or PETSC_DECIDE to have calculated)
79747c6ae99SBarry Smith .  dof - number of degrees of freedom per node
79847c6ae99SBarry Smith .  s - stencil width
79947c6ae99SBarry Smith -  lx, ly - arrays containing the number of nodes in each cell along
8000298fd71SBarry Smith            the x and y coordinates, or NULL. If non-null, these
80147c6ae99SBarry Smith            must be of length as m and n, and the corresponding
80247c6ae99SBarry Smith            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
80347c6ae99SBarry Smith            must be M, and the sum of the ly[] entries must be N.
80447c6ae99SBarry Smith 
80547c6ae99SBarry Smith    Output Parameter:
80647c6ae99SBarry Smith .  da - the resulting distributed array object
80747c6ae99SBarry Smith 
80847c6ae99SBarry Smith    Options Database Key:
809706a11cbSBarry Smith +  -dm_view - Calls DMView() at the conclusion of DMDACreate2d()
81047c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
81147c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
81247c6ae99SBarry Smith .  -da_processors_x <nx> - number of processors in x direction
81347c6ae99SBarry Smith .  -da_processors_y <ny> - number of processors in y direction
814e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement ratio in x direction
815e0f5d30fSBarry Smith .  -da_refine_y <ry> - refinement ratio in y direction
816e0f5d30fSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating, if M or N < 0
817e0f5d30fSBarry Smith 
81847c6ae99SBarry Smith 
81947c6ae99SBarry Smith    Level: beginner
82047c6ae99SBarry Smith 
82147c6ae99SBarry Smith    Notes:
822aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
823aa219208SBarry Smith    standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
82447c6ae99SBarry Smith    the standard 9-pt stencil.
82547c6ae99SBarry Smith 
826aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
827564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
828564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
82947c6ae99SBarry Smith 
830897f7067SBarry Smith    You must call DMSetUp() after this call before using this DM.
831897f7067SBarry Smith 
832897f7067SBarry Smith    If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call
833897f7067SBarry Smith    but before DMSetUp().
834897f7067SBarry Smith 
83547c6ae99SBarry Smith .keywords: distributed array, create, two-dimensional
83647c6ae99SBarry Smith 
837aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
83899f0b487SRichard Tran Mills           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(),
839d461ba97SBarry Smith           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
84047c6ae99SBarry Smith 
84147c6ae99SBarry Smith @*/
842fe16a2e9SBarry Smith 
843bff4a2f0SMatthew G. Knepley PetscErrorCode  DMDACreate2d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMDAStencilType stencil_type,
8449a42bb27SBarry Smith                           PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
84547c6ae99SBarry Smith {
84647c6ae99SBarry Smith   PetscErrorCode ierr;
84747c6ae99SBarry Smith 
84847c6ae99SBarry Smith   PetscFunctionBegin;
849aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
850c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*da, 2);CHKERRQ(ierr);
851aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, 1);CHKERRQ(ierr);
852aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr);
853bff4a2f0SMatthew G. Knepley   ierr = DMDASetBoundaryType(*da, bx, by, DM_BOUNDARY_NONE);CHKERRQ(ierr);
854aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
855aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
856aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
8570298fd71SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, NULL);CHKERRQ(ierr);
85847c6ae99SBarry Smith   PetscFunctionReturn(0);
85947c6ae99SBarry Smith }
860