xref: /petsc/src/dm/impls/da/da2.c (revision 7a8be3513cf479fb6a672bd91de7eb6883fcfd02)
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;
9c9493c35SStefano Zampini   PetscBool      iascii,isdraw,isglvis,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;
16ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRMPI(ierr);
1747c6ae99SBarry Smith 
18251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
19251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
20c9493c35SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr);
21251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
229a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
23251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr);
249a42bb27SBarry Smith #endif
2547c6ae99SBarry Smith   if (iascii) {
2647c6ae99SBarry Smith     PetscViewerFormat format;
2747c6ae99SBarry Smith 
2847c6ae99SBarry Smith     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
2976a8abe0SBarry Smith     if (format == PETSC_VIEWER_LOAD_BALANCE) {
3076a8abe0SBarry Smith       PetscInt      i,nmax = 0,nmin = PETSC_MAX_INT,navg = 0,*nz,nzlocal;
3176a8abe0SBarry Smith       DMDALocalInfo info;
3276a8abe0SBarry Smith       PetscMPIInt   size;
33ffc4695bSBarry Smith       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRMPI(ierr);
3476a8abe0SBarry Smith       ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
3576a8abe0SBarry Smith       nzlocal = info.xm*info.ym;
3676a8abe0SBarry Smith       ierr = PetscMalloc1(size,&nz);CHKERRQ(ierr);
37ffc4695bSBarry Smith       ierr = MPI_Allgather(&nzlocal,1,MPIU_INT,nz,1,MPIU_INT,PetscObjectComm((PetscObject)da));CHKERRMPI(ierr);
3876a8abe0SBarry Smith       for (i=0; i<(PetscInt)size; i++) {
3976a8abe0SBarry Smith         nmax = PetscMax(nmax,nz[i]);
4076a8abe0SBarry Smith         nmin = PetscMin(nmin,nz[i]);
4176a8abe0SBarry Smith         navg += nz[i];
4276a8abe0SBarry Smith       }
4376a8abe0SBarry Smith       ierr = PetscFree(nz);CHKERRQ(ierr);
4476a8abe0SBarry Smith       navg = navg/size;
4576a8abe0SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Load Balance - Grid Points: Min %D  avg %D  max %D\n",nmin,navg,nmax);CHKERRQ(ierr);
4676a8abe0SBarry Smith       PetscFunctionReturn(0);
4776a8abe0SBarry Smith     }
488ec8862eSJed Brown     if (format != PETSC_VIEWER_ASCII_VTK_DEPRECATED && format != PETSC_VIEWER_ASCII_VTK_CELL_DEPRECATED && format != PETSC_VIEWER_ASCII_GLVIS) {
49aa219208SBarry Smith       DMDALocalInfo info;
50aa219208SBarry Smith       ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
511575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
5247c6ae99SBarry 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);
5347c6ae99SBarry 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);
5447c6ae99SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
551575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
568135c375SStefano Zampini     } else if (format == PETSC_VIEWER_ASCII_GLVIS) {
578135c375SStefano Zampini       ierr = DMView_DA_GLVis(da,viewer);CHKERRQ(ierr);
583da9ae13SJed Brown     } else {
593da9ae13SJed Brown       ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr);
6047c6ae99SBarry Smith     }
6147c6ae99SBarry Smith   } else if (isdraw) {
6247c6ae99SBarry Smith     PetscDraw      draw;
6347c6ae99SBarry Smith     double         ymin = -1*dd->s-1,ymax = dd->N+dd->s;
6447c6ae99SBarry Smith     double         xmin = -1*dd->s-1,xmax = dd->M+dd->s;
6547c6ae99SBarry Smith     double         x,y;
668ea3bf28SBarry Smith     PetscInt       base;
678ea3bf28SBarry Smith     const PetscInt *idx;
6847c6ae99SBarry Smith     char           node[10];
6947c6ae99SBarry Smith     PetscBool      isnull;
7047c6ae99SBarry Smith 
7147c6ae99SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
725b399a63SLisandro Dalcin     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
735b399a63SLisandro Dalcin     if (isnull) PetscFunctionReturn(0);
7447c6ae99SBarry Smith 
755b399a63SLisandro Dalcin     ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
765b399a63SLisandro Dalcin     ierr = PetscDrawClear(draw);CHKERRQ(ierr);
775b399a63SLisandro Dalcin     ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
785b399a63SLisandro Dalcin 
795b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
8047c6ae99SBarry Smith     /* first processor draw all node lines */
81dd400576SPatrick Sanan     if (rank == 0) {
8247c6ae99SBarry Smith       ymin = 0.0; ymax = dd->N - 1;
8347c6ae99SBarry Smith       for (xmin=0; xmin<dd->M; xmin++) {
8447c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
8547c6ae99SBarry Smith       }
8647c6ae99SBarry Smith       xmin = 0.0; xmax = dd->M - 1;
8747c6ae99SBarry Smith       for (ymin=0; ymin<dd->N; ymin++) {
8847c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
8947c6ae99SBarry Smith       }
9047c6ae99SBarry Smith     }
915b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
925b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
9347c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
9447c6ae99SBarry Smith 
955b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
9647c6ae99SBarry Smith     /* draw my box */
975b399a63SLisandro Dalcin     xmin = dd->xs/dd->w; xmax =(dd->xe-1)/dd->w; ymin = dd->ys; ymax = dd->ye - 1;
9847c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
9947c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
10047c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
10147c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
10247c6ae99SBarry Smith     /* put in numbers */
10347c6ae99SBarry Smith     base = (dd->base)/dd->w;
10447c6ae99SBarry Smith     for (y=ymin; y<=ymax; y++) {
10547c6ae99SBarry Smith       for (x=xmin; x<=xmax; x++) {
1065b399a63SLisandro Dalcin         ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)base++);CHKERRQ(ierr);
10747c6ae99SBarry Smith         ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
10847c6ae99SBarry Smith       }
10947c6ae99SBarry Smith     }
1105b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
1115b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
11247c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
11347c6ae99SBarry Smith 
1145b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
1155b399a63SLisandro Dalcin     /* overlay ghost numbers, useful for error checking */
11645b6f7e9SBarry Smith     ierr = ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr);
1175b399a63SLisandro Dalcin     base = 0; xmin = dd->Xs; xmax = dd->Xe; ymin = dd->Ys; ymax = dd->Ye;
11847c6ae99SBarry Smith     for (y=ymin; y<ymax; y++) {
11947c6ae99SBarry Smith       for (x=xmin; x<xmax; x++) {
12047c6ae99SBarry Smith         if ((base % dd->w) == 0) {
1215b399a63SLisandro Dalcin           ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)(idx[base/dd->w]));CHKERRQ(ierr);
12247c6ae99SBarry Smith           ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
12347c6ae99SBarry Smith         }
12447c6ae99SBarry Smith         base++;
12547c6ae99SBarry Smith       }
12647c6ae99SBarry Smith     }
127302440fdSBarry Smith     ierr = ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr);
1285b399a63SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
1295b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
13047c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
131832b7cebSLisandro Dalcin     ierr = PetscDrawSave(draw);CHKERRQ(ierr);
132c9493c35SStefano Zampini   } else if (isglvis) {
133c9493c35SStefano Zampini     ierr = DMView_DA_GLVis(da,viewer);CHKERRQ(ierr);
1349a42bb27SBarry Smith   } else if (isbinary) {
1359a42bb27SBarry Smith     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
1369a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1379a42bb27SBarry Smith   } else if (ismatlab) {
1389a42bb27SBarry Smith     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
1399a42bb27SBarry Smith #endif
14011aeaf0aSBarry Smith   }
14147c6ae99SBarry Smith   PetscFunctionReturn(0);
14247c6ae99SBarry Smith }
14347c6ae99SBarry Smith 
14447c6ae99SBarry Smith #if defined(new)
14547c6ae99SBarry Smith /*
146aa219208SBarry Smith   DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
147aa219208SBarry Smith     function lives on a DMDA
14847c6ae99SBarry Smith 
14947c6ae99SBarry Smith         y ~= (F(u + ha) - F(u))/h,
15047c6ae99SBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
15147c6ae99SBarry Smith         u = current iterate
15247c6ae99SBarry Smith         h = difference interval
15347c6ae99SBarry Smith */
154aa219208SBarry Smith PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a)
15547c6ae99SBarry Smith {
15647c6ae99SBarry Smith   PetscScalar    h,*aa,*ww,v;
15747c6ae99SBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
15847c6ae99SBarry Smith   PetscErrorCode ierr;
15947c6ae99SBarry Smith   PetscInt       gI,nI;
16047c6ae99SBarry Smith   MatStencil     stencil;
161aa219208SBarry Smith   DMDALocalInfo  info;
16247c6ae99SBarry Smith 
16347c6ae99SBarry Smith   PetscFunctionBegin;
16447c6ae99SBarry Smith   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
16547c6ae99SBarry Smith   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
16647c6ae99SBarry Smith 
16747c6ae99SBarry Smith   ierr = VecGetArray(U,&ww);CHKERRQ(ierr);
16847c6ae99SBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
16947c6ae99SBarry Smith 
17047c6ae99SBarry Smith   nI = 0;
17147c6ae99SBarry Smith   h  = ww[gI];
17247c6ae99SBarry Smith   if (h == 0.0) h = 1.0;
17347c6ae99SBarry Smith   if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin;
17447c6ae99SBarry Smith   else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
17547c6ae99SBarry Smith   h *= epsilon;
17647c6ae99SBarry Smith 
17747c6ae99SBarry Smith   ww[gI] += h;
17847c6ae99SBarry Smith   ierr    = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
17947c6ae99SBarry Smith   aa[nI]  = (v - aa[nI])/h;
18047c6ae99SBarry Smith   ww[gI] -= h;
18147c6ae99SBarry Smith   nI++;
1828865f1eaSKarl Rupp 
18347c6ae99SBarry Smith   ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr);
18447c6ae99SBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
18547c6ae99SBarry Smith   PetscFunctionReturn(0);
18647c6ae99SBarry Smith }
18747c6ae99SBarry Smith #endif
18847c6ae99SBarry Smith 
1897087cfbeSBarry Smith PetscErrorCode  DMSetUp_DA_2D(DM da)
19047c6ae99SBarry Smith {
19147c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
19247c6ae99SBarry Smith   const PetscInt   M            = dd->M;
19347c6ae99SBarry Smith   const PetscInt   N            = dd->N;
19447c6ae99SBarry Smith   PetscInt         m            = dd->m;
19547c6ae99SBarry Smith   PetscInt         n            = dd->n;
19647c6ae99SBarry Smith   const PetscInt   dof          = dd->w;
19747c6ae99SBarry Smith   const PetscInt   s            = dd->s;
198bff4a2f0SMatthew G. Knepley   DMBoundaryType   bx           = dd->bx;
199bff4a2f0SMatthew G. Knepley   DMBoundaryType   by           = dd->by;
20019fd82e9SBarry Smith   DMDAStencilType  stencil_type = dd->stencil_type;
20147c6ae99SBarry Smith   PetscInt         *lx          = dd->lx;
20247c6ae99SBarry Smith   PetscInt         *ly          = dd->ly;
20347c6ae99SBarry Smith   MPI_Comm         comm;
20447c6ae99SBarry Smith   PetscMPIInt      rank,size;
205bd1fc5aeSBarry Smith   PetscInt         xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,IXs,IXe,IYs,IYe;
2068ea3bf28SBarry Smith   PetscInt         up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn;
207db87c5ecSEthan Coon   PetscInt         xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
20847c6ae99SBarry Smith   PetscInt         s_x,s_y; /* s proportionalized to w */
20947c6ae99SBarry Smith   PetscInt         sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
21047c6ae99SBarry Smith   Vec              local,global;
211bd1fc5aeSBarry Smith   VecScatter       gtol;
21245b6f7e9SBarry Smith   IS               to,from;
21347c6ae99SBarry Smith   PetscErrorCode   ierr;
21447c6ae99SBarry Smith 
21547c6ae99SBarry Smith   PetscFunctionBegin;
216*7a8be351SBarry Smith   PetscCheck(stencil_type != DMDA_STENCIL_BOX || (bx != DM_BOUNDARY_MIRROR && by != DM_BOUNDARY_MIRROR),PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Mirror boundary and box stencil");
21747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
2183855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES)
219*7a8be351SBarry Smith   PetscCheck(((PetscInt64) M)*((PetscInt64) N)*((PetscInt64) dof) <= (PetscInt64) PETSC_MPI_INT_MAX,comm,PETSC_ERR_INT_OVERFLOW,"Mesh of %D by %D by %D (dof) is too large for 32 bit indices",M,N,dof);
2203855c12bSBarry Smith #endif
22147c6ae99SBarry Smith 
222ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
223ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
22447c6ae99SBarry Smith 
2257d310018SBarry Smith   dd->p = 1;
22647c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
227*7a8be351SBarry Smith     PetscCheck(m >= 1,comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
228*7a8be351SBarry Smith     else PetscCheck(m <= size,comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
22947c6ae99SBarry Smith   }
23047c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
231*7a8be351SBarry Smith     PetscCheck(n >= 1,comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
232*7a8be351SBarry Smith     else PetscCheck(n <= size,comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
23347c6ae99SBarry Smith   }
23447c6ae99SBarry Smith 
23547c6ae99SBarry Smith   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
23647c6ae99SBarry Smith     if (n != PETSC_DECIDE) {
23747c6ae99SBarry Smith       m = size/n;
23847c6ae99SBarry Smith     } else if (m != PETSC_DECIDE) {
23947c6ae99SBarry Smith       n = size/m;
24047c6ae99SBarry Smith     } else {
24147c6ae99SBarry Smith       /* try for squarish distribution */
242369cc0aeSBarry Smith       m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N)));
24347c6ae99SBarry Smith       if (!m) m = 1;
24447c6ae99SBarry Smith       while (m > 0) {
24547c6ae99SBarry Smith         n = size/m;
24647c6ae99SBarry Smith         if (m*n == size) break;
24747c6ae99SBarry Smith         m--;
24847c6ae99SBarry Smith       }
24947c6ae99SBarry Smith       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
25047c6ae99SBarry Smith     }
251*7a8be351SBarry Smith     PetscCheck(m*n == size,comm,PETSC_ERR_PLIB,"Unable to create partition, check the size of the communicator and input m and n ");
252*7a8be351SBarry Smith   } else PetscCheck(m*n == size,comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
25347c6ae99SBarry Smith 
254*7a8be351SBarry Smith   PetscCheck(M >= m,comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
255*7a8be351SBarry Smith   PetscCheck(N >= n,comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
25647c6ae99SBarry Smith 
25747c6ae99SBarry Smith   /*
25847c6ae99SBarry Smith      Determine locally owned region
25947c6ae99SBarry Smith      xs is the first local node number, x is the number of local nodes
26047c6ae99SBarry Smith   */
26147c6ae99SBarry Smith   if (!lx) {
262785e854fSJed Brown     ierr = PetscMalloc1(m, &dd->lx);CHKERRQ(ierr);
26347c6ae99SBarry Smith     lx   = dd->lx;
26447c6ae99SBarry Smith     for (i=0; i<m; i++) {
26547c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > i);
26647c6ae99SBarry Smith     }
26747c6ae99SBarry Smith   }
26847c6ae99SBarry Smith   x  = lx[rank % m];
26947c6ae99SBarry Smith   xs = 0;
27047c6ae99SBarry Smith   for (i=0; i<(rank % m); i++) {
27147c6ae99SBarry Smith     xs += lx[i];
27247c6ae99SBarry Smith   }
27376bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
27447c6ae99SBarry Smith     left = xs;
27547c6ae99SBarry Smith     for (i=(rank % m); i<m; i++) {
27647c6ae99SBarry Smith       left += lx[i];
27747c6ae99SBarry Smith     }
278*7a8be351SBarry Smith     PetscCheck(left == M,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M: %D %D",left,M);
27976bd3646SJed Brown   }
28047c6ae99SBarry Smith 
28147c6ae99SBarry Smith   /*
28247c6ae99SBarry Smith      Determine locally owned region
28347c6ae99SBarry Smith      ys is the first local node number, y is the number of local nodes
28447c6ae99SBarry Smith   */
28547c6ae99SBarry Smith   if (!ly) {
286785e854fSJed Brown     ierr = PetscMalloc1(n, &dd->ly);CHKERRQ(ierr);
28747c6ae99SBarry Smith     ly   = dd->ly;
28847c6ae99SBarry Smith     for (i=0; i<n; i++) {
28947c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > i);
29047c6ae99SBarry Smith     }
29147c6ae99SBarry Smith   }
29247c6ae99SBarry Smith   y  = ly[rank/m];
29347c6ae99SBarry Smith   ys = 0;
29447c6ae99SBarry Smith   for (i=0; i<(rank/m); i++) {
29547c6ae99SBarry Smith     ys += ly[i];
29647c6ae99SBarry Smith   }
29776bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
29847c6ae99SBarry Smith     left = ys;
29947c6ae99SBarry Smith     for (i=(rank/m); i<n; i++) {
30047c6ae99SBarry Smith       left += ly[i];
30147c6ae99SBarry Smith     }
302*7a8be351SBarry Smith     PetscCheck(left == N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of ly across processors not equal to N: %D %D",left,N);
30376bd3646SJed Brown   }
30447c6ae99SBarry Smith 
305bcea557cSEthan Coon   /*
306bcea557cSEthan Coon    check if the scatter requires more than one process neighbor or wraps around
307bcea557cSEthan Coon    the domain more than once
308bcea557cSEthan Coon   */
309*7a8be351SBarry Smith   PetscCheck((x >= s) || ((m <= 1) && (bx != DM_BOUNDARY_PERIODIC)),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s);
310*7a8be351SBarry Smith   PetscCheck((y >= s) || ((n <= 1) && (by != DM_BOUNDARY_PERIODIC)),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s);
31147c6ae99SBarry Smith   xe = xs + x;
31247c6ae99SBarry Smith   ye = ys + y;
31347c6ae99SBarry Smith 
314ce00eea3SSatish Balay   /* determine ghost region (Xs) and region scattered into (IXs)  */
315d9c9ebe5SPeter Brune   if (xs-s > 0) {
316d9c9ebe5SPeter Brune     Xs = xs - s; IXs = xs - s;
31788661749SPeter Brune   } else {
31888661749SPeter Brune     if (bx) {
31988661749SPeter Brune       Xs = xs - s;
32088661749SPeter Brune     } else {
32188661749SPeter Brune       Xs = 0;
32288661749SPeter Brune     }
32388661749SPeter Brune     IXs = 0;
32488661749SPeter Brune   }
325d9c9ebe5SPeter Brune   if (xe+s <= M) {
326d9c9ebe5SPeter Brune     Xe = xe + s; IXe = xe + s;
32788661749SPeter Brune   } else {
32888661749SPeter Brune     if (bx) {
329d9c9ebe5SPeter Brune       Xs = xs - s; Xe = xe + s;
33088661749SPeter Brune     } else {
33188661749SPeter Brune       Xe = M;
33288661749SPeter Brune     }
33388661749SPeter Brune     IXe = M;
33488661749SPeter Brune   }
33547c6ae99SBarry Smith 
336bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) {
337d9c9ebe5SPeter Brune     IXs = xs - s;
338d9c9ebe5SPeter Brune     IXe = xe + s;
339d9c9ebe5SPeter Brune     Xs  = xs - s;
340d9c9ebe5SPeter Brune     Xe  = xe + s;
34188661749SPeter Brune   }
34247c6ae99SBarry Smith 
343d9c9ebe5SPeter Brune   if (ys-s > 0) {
344d9c9ebe5SPeter Brune     Ys = ys - s; IYs = ys - s;
34588661749SPeter Brune   } else {
34688661749SPeter Brune     if (by) {
34788661749SPeter Brune       Ys = ys - s;
34888661749SPeter Brune     } else {
34988661749SPeter Brune       Ys = 0;
35088661749SPeter Brune     }
35188661749SPeter Brune     IYs = 0;
35288661749SPeter Brune   }
353d9c9ebe5SPeter Brune   if (ye+s <= N) {
354d9c9ebe5SPeter Brune     Ye = ye + s; IYe = ye + s;
35588661749SPeter Brune   } else {
35688661749SPeter Brune     if (by) {
35788661749SPeter Brune       Ye = ye + s;
35888661749SPeter Brune     } else {
35988661749SPeter Brune       Ye = N;
36088661749SPeter Brune     }
36188661749SPeter Brune     IYe = N;
36288661749SPeter Brune   }
36388661749SPeter Brune 
364bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) {
365d9c9ebe5SPeter Brune     IYs = ys - s;
366d9c9ebe5SPeter Brune     IYe = ye + s;
367d9c9ebe5SPeter Brune     Ys  = ys - s;
368d9c9ebe5SPeter Brune     Ye  = ye + s;
36988661749SPeter Brune   }
37088661749SPeter Brune 
37188661749SPeter Brune   /* stencil length in each direction */
372d9c9ebe5SPeter Brune   s_x = s;
373d9c9ebe5SPeter Brune   s_y = s;
37447c6ae99SBarry Smith 
37547c6ae99SBarry Smith   /* determine starting point of each processor */
37647c6ae99SBarry Smith   nn       = x*y;
377dcca6d9dSJed Brown   ierr     = PetscMalloc2(size+1,&bases,size,&ldims);CHKERRQ(ierr);
37855b25c41SPierre Jolivet   ierr     = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRMPI(ierr);
37947c6ae99SBarry Smith   bases[0] = 0;
38047c6ae99SBarry Smith   for (i=1; i<=size; i++) {
38147c6ae99SBarry Smith     bases[i] = ldims[i-1];
38247c6ae99SBarry Smith   }
38347c6ae99SBarry Smith   for (i=1; i<=size; i++) {
38447c6ae99SBarry Smith     bases[i] += bases[i-1];
38547c6ae99SBarry Smith   }
386ce00eea3SSatish Balay   base = bases[rank]*dof;
38747c6ae99SBarry Smith 
38847c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
389ce00eea3SSatish Balay   dd->Nlocal = x*y*dof;
390b1fb7eb7SBarry Smith   ierr       = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global);CHKERRQ(ierr);
391ce00eea3SSatish Balay   dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof;
392b1fb7eb7SBarry Smith   ierr       = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local);CHKERRQ(ierr);
39347c6ae99SBarry Smith 
3944104a7a0SPatrick Sanan   /* generate global to local vector scatter and local to global mapping*/
39547c6ae99SBarry Smith 
396ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
397ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
3984104a7a0SPatrick Sanan   ierr  = PetscMalloc1((IXe-IXs)*(IYe-IYs),&idx);CHKERRQ(ierr);
399d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_BOX) {
400ce00eea3SSatish Balay     left  = IXs - Xs; right = left + (IXe-IXs);
401ce00eea3SSatish Balay     down  = IYs - Ys; up = down + (IYe-IYs);
402ce00eea3SSatish Balay     count = 0;
403ce00eea3SSatish Balay     for (i=down; i<up; i++) {
404ce00eea3SSatish Balay       for (j=left; j<right; j++) {
405ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
406ce00eea3SSatish Balay       }
407ce00eea3SSatish Balay     }
408ce00eea3SSatish Balay     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
409ce00eea3SSatish Balay 
41047c6ae99SBarry Smith   } else {
41147c6ae99SBarry Smith     /* must drop into cross shape region */
41247c6ae99SBarry Smith     /*       ---------|
41347c6ae99SBarry Smith             |  top    |
414ce00eea3SSatish Balay          |---         ---| up
41547c6ae99SBarry Smith          |   middle      |
41647c6ae99SBarry Smith          |               |
417ce00eea3SSatish Balay          ----         ---- down
41847c6ae99SBarry Smith             | bottom  |
41947c6ae99SBarry Smith             -----------
42047c6ae99SBarry Smith          Xs xs        xe Xe */
421ce00eea3SSatish Balay     left  = xs - Xs; right = left + x;
422ce00eea3SSatish Balay     down  = ys - Ys; up = down + y;
42347c6ae99SBarry Smith     count = 0;
424ce00eea3SSatish Balay     /* bottom */
425ce00eea3SSatish Balay     for (i=(IYs-Ys); i<down; i++) {
426ce00eea3SSatish Balay       for (j=left; j<right; j++) {
427ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
42847c6ae99SBarry Smith       }
42947c6ae99SBarry Smith     }
43047c6ae99SBarry Smith     /* middle */
43147c6ae99SBarry Smith     for (i=down; i<up; i++) {
432ce00eea3SSatish Balay       for (j=(IXs-Xs); j<(IXe-Xs); j++) {
433ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
43447c6ae99SBarry Smith       }
43547c6ae99SBarry Smith     }
43647c6ae99SBarry Smith     /* top */
437ce00eea3SSatish Balay     for (i=up; i<up+IYe-ye; i++) {
438ce00eea3SSatish Balay       for (j=left; j<right; j++) {
439ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
44047c6ae99SBarry Smith       }
44147c6ae99SBarry Smith     }
44247c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
44347c6ae99SBarry Smith   }
44447c6ae99SBarry Smith 
44547c6ae99SBarry Smith   /* determine who lies on each side of us stored in    n6 n7 n8
44647c6ae99SBarry Smith                                                         n3    n5
44747c6ae99SBarry Smith                                                         n0 n1 n2
44847c6ae99SBarry Smith   */
44947c6ae99SBarry Smith 
45047c6ae99SBarry Smith   /* Assume the Non-Periodic Case */
45147c6ae99SBarry Smith   n1 = rank - m;
45247c6ae99SBarry Smith   if (rank % m) {
45347c6ae99SBarry Smith     n0 = n1 - 1;
45447c6ae99SBarry Smith   } else {
45547c6ae99SBarry Smith     n0 = -1;
45647c6ae99SBarry Smith   }
45747c6ae99SBarry Smith   if ((rank+1) % m) {
45847c6ae99SBarry Smith     n2 = n1 + 1;
45947c6ae99SBarry Smith     n5 = rank + 1;
46047c6ae99SBarry Smith     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
46147c6ae99SBarry Smith   } else {
46247c6ae99SBarry Smith     n2 = -1; n5 = -1; n8 = -1;
46347c6ae99SBarry Smith   }
46447c6ae99SBarry Smith   if (rank % m) {
46547c6ae99SBarry Smith     n3 = rank - 1;
46647c6ae99SBarry Smith     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
46747c6ae99SBarry Smith   } else {
46847c6ae99SBarry Smith     n3 = -1; n6 = -1;
46947c6ae99SBarry Smith   }
47047c6ae99SBarry Smith   n7 = rank + m; if (n7 >= m*n) n7 = -1;
47147c6ae99SBarry Smith 
472bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC && by == DM_BOUNDARY_PERIODIC) {
47347c6ae99SBarry Smith     /* Modify for Periodic Cases */
47447c6ae99SBarry Smith     /* Handle all four corners */
47547c6ae99SBarry Smith     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
47647c6ae99SBarry Smith     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
47747c6ae99SBarry Smith     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
47847c6ae99SBarry Smith     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
47947c6ae99SBarry Smith 
48047c6ae99SBarry Smith     /* Handle Top and Bottom Sides */
48147c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
48247c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
48347c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
48447c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
48547c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
48647c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
48747c6ae99SBarry Smith 
48847c6ae99SBarry Smith     /* Handle Left and Right Sides */
48947c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
49047c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
49147c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
49247c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
49347c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
49447c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
495bff4a2f0SMatthew G. Knepley   } else if (by == DM_BOUNDARY_PERIODIC) {  /* Handle Top and Bottom Sides */
496ce00eea3SSatish Balay     if (n1 < 0) n1 = rank + m * (n-1);
497ce00eea3SSatish Balay     if (n7 < 0) n7 = rank - m * (n-1);
498ce00eea3SSatish Balay     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
499ce00eea3SSatish Balay     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
500ce00eea3SSatish Balay     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
501ce00eea3SSatish Balay     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
502bff4a2f0SMatthew G. Knepley   } else if (bx == DM_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */
503ce00eea3SSatish Balay     if (n3 < 0) n3 = rank + (m-1);
504ce00eea3SSatish Balay     if (n5 < 0) n5 = rank - (m-1);
505ce00eea3SSatish Balay     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
506ce00eea3SSatish Balay     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
507ce00eea3SSatish Balay     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
508ce00eea3SSatish Balay     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
50947c6ae99SBarry Smith   }
510ce00eea3SSatish Balay 
511785e854fSJed Brown   ierr = PetscMalloc1(9,&dd->neighbors);CHKERRQ(ierr);
5128865f1eaSKarl Rupp 
51347c6ae99SBarry Smith   dd->neighbors[0] = n0;
51447c6ae99SBarry Smith   dd->neighbors[1] = n1;
51547c6ae99SBarry Smith   dd->neighbors[2] = n2;
51647c6ae99SBarry Smith   dd->neighbors[3] = n3;
51747c6ae99SBarry Smith   dd->neighbors[4] = rank;
51847c6ae99SBarry Smith   dd->neighbors[5] = n5;
51947c6ae99SBarry Smith   dd->neighbors[6] = n6;
52047c6ae99SBarry Smith   dd->neighbors[7] = n7;
52147c6ae99SBarry Smith   dd->neighbors[8] = n8;
52247c6ae99SBarry Smith 
523d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
52447c6ae99SBarry Smith     /* save corner processor numbers */
52547c6ae99SBarry Smith     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
52647c6ae99SBarry Smith     n0  = n2 = n6 = n8 = -1;
52747c6ae99SBarry Smith   }
52847c6ae99SBarry Smith 
529785e854fSJed Brown   ierr = PetscMalloc1((Xe-Xs)*(Ye-Ys),&idx);CHKERRQ(ierr);
53047c6ae99SBarry Smith 
531ce00eea3SSatish Balay   nn = 0;
53247c6ae99SBarry Smith   xbase = bases[rank];
53347c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
53447c6ae99SBarry Smith     if (n0 >= 0) { /* left below */
535ce00eea3SSatish Balay       x_t = lx[n0 % m];
53647c6ae99SBarry Smith       y_t = ly[(n0/m)];
53747c6ae99SBarry Smith       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
5388865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
53947c6ae99SBarry Smith     }
540ac119b13SBarry Smith 
54147c6ae99SBarry Smith     if (n1 >= 0) { /* directly below */
54247c6ae99SBarry Smith       x_t = x;
54347c6ae99SBarry Smith       y_t = ly[(n1/m)];
54447c6ae99SBarry Smith       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
5458865f1eaSKarl Rupp       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
546bff4a2f0SMatthew G. Knepley     } else if (by == DM_BOUNDARY_MIRROR) {
5478865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
54847c6ae99SBarry Smith     }
549ac119b13SBarry Smith 
55047c6ae99SBarry Smith     if (n2 >= 0) { /* right below */
551ce00eea3SSatish Balay       x_t = lx[n2 % m];
55247c6ae99SBarry Smith       y_t = ly[(n2/m)];
55347c6ae99SBarry Smith       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
5548865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
55547c6ae99SBarry Smith     }
55647c6ae99SBarry Smith   }
55747c6ae99SBarry Smith 
55847c6ae99SBarry Smith   for (i=0; i<y; i++) {
55947c6ae99SBarry Smith     if (n3 >= 0) { /* directly left */
560ce00eea3SSatish Balay       x_t = lx[n3 % m];
56147c6ae99SBarry Smith       /* y_t = y; */
56247c6ae99SBarry Smith       s_t = bases[n3] + (i+1)*x_t - s_x;
5638865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
564bff4a2f0SMatthew G. Knepley     } else if (bx == DM_BOUNDARY_MIRROR) {
5658865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
56647c6ae99SBarry Smith     }
56747c6ae99SBarry Smith 
5688865f1eaSKarl Rupp     for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
56947c6ae99SBarry Smith 
57047c6ae99SBarry Smith     if (n5 >= 0) { /* directly right */
571ce00eea3SSatish Balay       x_t = lx[n5 % m];
57247c6ae99SBarry Smith       /* y_t = y; */
57347c6ae99SBarry Smith       s_t = bases[n5] + (i)*x_t;
5748865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
575bff4a2f0SMatthew G. Knepley     } else if (bx == DM_BOUNDARY_MIRROR) {
5768865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
57747c6ae99SBarry Smith     }
57847c6ae99SBarry Smith   }
57947c6ae99SBarry Smith 
58047c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
58147c6ae99SBarry Smith     if (n6 >= 0) { /* left above */
582ce00eea3SSatish Balay       x_t = lx[n6 % m];
58347c6ae99SBarry Smith       /* y_t = ly[(n6/m)]; */
58447c6ae99SBarry Smith       s_t = bases[n6] + (i)*x_t - s_x;
5858865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
58647c6ae99SBarry Smith     }
587ac119b13SBarry Smith 
58847c6ae99SBarry Smith     if (n7 >= 0) { /* directly above */
58947c6ae99SBarry Smith       x_t = x;
59047c6ae99SBarry Smith       /* y_t = ly[(n7/m)]; */
59147c6ae99SBarry Smith       s_t = bases[n7] + (i-1)*x_t;
5928865f1eaSKarl Rupp       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
593bff4a2f0SMatthew G. Knepley     } else if (by == DM_BOUNDARY_MIRROR) {
5948865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
59547c6ae99SBarry Smith     }
596ac119b13SBarry Smith 
59747c6ae99SBarry Smith     if (n8 >= 0) { /* right above */
598ce00eea3SSatish Balay       x_t = lx[n8 % m];
59947c6ae99SBarry Smith       /* y_t = ly[(n8/m)]; */
60047c6ae99SBarry Smith       s_t = bases[n8] + (i-1)*x_t;
6018865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
60247c6ae99SBarry Smith     }
60347c6ae99SBarry Smith   }
60447c6ae99SBarry Smith 
605b1fb7eb7SBarry Smith   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from);CHKERRQ(ierr);
6069448b7f1SJunchao Zhang   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
6073bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr);
608fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
609fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
61047c6ae99SBarry Smith 
611d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
612ce00eea3SSatish Balay     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
613ce00eea3SSatish Balay   }
614ce00eea3SSatish Balay 
615288e7d53SBarry Smith   if (((stencil_type == DMDA_STENCIL_STAR)  || (bx && bx != DM_BOUNDARY_PERIODIC) || (by && by != DM_BOUNDARY_PERIODIC))) {
61647c6ae99SBarry Smith     /*
61747c6ae99SBarry Smith         Recompute the local to global mappings, this time keeping the
618ce00eea3SSatish Balay       information about the cross corner processor numbers and any ghosted
619ce00eea3SSatish Balay       but not periodic indices.
62047c6ae99SBarry Smith     */
62147c6ae99SBarry Smith     nn    = 0;
62247c6ae99SBarry Smith     xbase = bases[rank];
62347c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
62447c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
625ce00eea3SSatish Balay         x_t = lx[n0 % m];
62647c6ae99SBarry Smith         y_t = ly[(n0/m)];
62747c6ae99SBarry Smith         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
6288865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
629ce00eea3SSatish Balay       } else if (xs-Xs > 0 && ys-Ys > 0) {
6308865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
63147c6ae99SBarry Smith       }
63247c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
63347c6ae99SBarry Smith         x_t = x;
63447c6ae99SBarry Smith         y_t = ly[(n1/m)];
63547c6ae99SBarry Smith         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
6368865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
637ce00eea3SSatish Balay       } else if (ys-Ys > 0) {
638bff4a2f0SMatthew G. Knepley         if (by == DM_BOUNDARY_MIRROR) {
6398865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
640624904c4SBarry Smith         } else {
6418865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
64247c6ae99SBarry Smith         }
643624904c4SBarry Smith       }
64447c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
645ce00eea3SSatish Balay         x_t = lx[n2 % m];
64647c6ae99SBarry Smith         y_t = ly[(n2/m)];
64747c6ae99SBarry Smith         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
6488865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
649ce00eea3SSatish Balay       } else if (Xe-xe> 0 && ys-Ys > 0) {
6508865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
65147c6ae99SBarry Smith       }
65247c6ae99SBarry Smith     }
65347c6ae99SBarry Smith 
65447c6ae99SBarry Smith     for (i=0; i<y; i++) {
65547c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
656ce00eea3SSatish Balay         x_t = lx[n3 % m];
65747c6ae99SBarry Smith         /* y_t = y; */
65847c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x;
6598865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
660ce00eea3SSatish Balay       } else if (xs-Xs > 0) {
661bff4a2f0SMatthew G. Knepley         if (bx == DM_BOUNDARY_MIRROR) {
6628865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
663624904c4SBarry Smith         } else {
6648865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
66547c6ae99SBarry Smith         }
666624904c4SBarry Smith       }
66747c6ae99SBarry Smith 
6688865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
66947c6ae99SBarry Smith 
67047c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
671ce00eea3SSatish Balay         x_t = lx[n5 % m];
67247c6ae99SBarry Smith         /* y_t = y; */
67347c6ae99SBarry Smith         s_t = bases[n5] + (i)*x_t;
6748865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
675ce00eea3SSatish Balay       } else if (Xe-xe > 0) {
676bff4a2f0SMatthew G. Knepley         if (bx == DM_BOUNDARY_MIRROR) {
6778865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
678624904c4SBarry Smith         } else {
6798865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
68047c6ae99SBarry Smith         }
68147c6ae99SBarry Smith       }
682624904c4SBarry Smith     }
68347c6ae99SBarry Smith 
68447c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
68547c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
686ce00eea3SSatish Balay         x_t = lx[n6 % m];
68747c6ae99SBarry Smith         /* y_t = ly[(n6/m)]; */
68847c6ae99SBarry Smith         s_t = bases[n6] + (i)*x_t - s_x;
6898865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
690ce00eea3SSatish Balay       } else if (xs-Xs > 0 && Ye-ye > 0) {
6918865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
69247c6ae99SBarry Smith       }
69347c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
69447c6ae99SBarry Smith         x_t = x;
69547c6ae99SBarry Smith         /* y_t = ly[(n7/m)]; */
69647c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t;
6978865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
698ce00eea3SSatish Balay       } else if (Ye-ye > 0) {
699bff4a2f0SMatthew G. Knepley         if (by == DM_BOUNDARY_MIRROR) {
7008865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
701624904c4SBarry Smith         } else {
7028865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
70347c6ae99SBarry Smith         }
704624904c4SBarry Smith       }
70547c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
706ce00eea3SSatish Balay         x_t = lx[n8 % m];
70747c6ae99SBarry Smith         /* y_t = ly[(n8/m)]; */
70847c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t;
7098865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
710ce00eea3SSatish Balay       } else if (Xe-xe > 0 && Ye-ye > 0) {
7118865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
71247c6ae99SBarry Smith       }
71347c6ae99SBarry Smith     }
71447c6ae99SBarry Smith   }
715ce00eea3SSatish Balay   /*
716ce00eea3SSatish Balay      Set the local to global ordering in the global vector, this allows use
717ce00eea3SSatish Balay      of VecSetValuesLocal().
718ce00eea3SSatish Balay   */
71945b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr);
7203bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr);
72147c6ae99SBarry Smith 
722ce00eea3SSatish Balay   ierr  = PetscFree2(bases,ldims);CHKERRQ(ierr);
72347c6ae99SBarry Smith   dd->m = m;  dd->n  = n;
724ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
725ce00eea3SSatish Balay   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
726ce00eea3SSatish Balay   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
72747c6ae99SBarry Smith 
728fcfd50ebSBarry Smith   ierr = VecDestroy(&local);CHKERRQ(ierr);
729fcfd50ebSBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
73047c6ae99SBarry Smith 
73147c6ae99SBarry Smith   dd->gtol      = gtol;
73247c6ae99SBarry Smith   dd->base      = base;
7339a42bb27SBarry Smith   da->ops->view = DMView_DA_2d;
7340298fd71SBarry Smith   dd->ltol      = NULL;
7350298fd71SBarry Smith   dd->ao        = NULL;
73647c6ae99SBarry Smith   PetscFunctionReturn(0);
73747c6ae99SBarry Smith }
73847c6ae99SBarry Smith 
73947c6ae99SBarry Smith /*@C
740aa219208SBarry Smith    DMDACreate2d -  Creates an object that will manage the communication of  two-dimensional
74147c6ae99SBarry Smith    regular array data that is distributed across some processors.
74247c6ae99SBarry Smith 
743d083f849SBarry Smith    Collective
74447c6ae99SBarry Smith 
74547c6ae99SBarry Smith    Input Parameters:
74647c6ae99SBarry Smith +  comm - MPI communicator
7471321219cSEthan Coon .  bx,by - type of ghost nodes the array have.
748bff4a2f0SMatthew G. Knepley          Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC.
749aa219208SBarry Smith .  stencil_type - stencil type.  Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
750897f7067SBarry Smith .  M,N - global dimension in each direction of the array
75147c6ae99SBarry Smith .  m,n - corresponding number of processors in each dimension
75247c6ae99SBarry Smith          (or PETSC_DECIDE to have calculated)
75347c6ae99SBarry Smith .  dof - number of degrees of freedom per node
75447c6ae99SBarry Smith .  s - stencil width
75547c6ae99SBarry Smith -  lx, ly - arrays containing the number of nodes in each cell along
7560298fd71SBarry Smith            the x and y coordinates, or NULL. If non-null, these
75747c6ae99SBarry Smith            must be of length as m and n, and the corresponding
75847c6ae99SBarry Smith            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
75947c6ae99SBarry Smith            must be M, and the sum of the ly[] entries must be N.
76047c6ae99SBarry Smith 
76147c6ae99SBarry Smith    Output Parameter:
76247c6ae99SBarry Smith .  da - the resulting distributed array object
76347c6ae99SBarry Smith 
76447c6ae99SBarry Smith    Options Database Key:
765706a11cbSBarry Smith +  -dm_view - Calls DMView() at the conclusion of DMDACreate2d()
766e43dc8daSBarry Smith .  -da_grid_x <nx> - number of grid points in x direction
767e43dc8daSBarry Smith .  -da_grid_y <ny> - number of grid points in y direction
76847c6ae99SBarry Smith .  -da_processors_x <nx> - number of processors in x direction
76947c6ae99SBarry Smith .  -da_processors_y <ny> - number of processors in y direction
770e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement ratio in x direction
771e0f5d30fSBarry Smith .  -da_refine_y <ry> - refinement ratio in y direction
772e43dc8daSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating
773e0f5d30fSBarry Smith 
77447c6ae99SBarry Smith    Level: beginner
77547c6ae99SBarry Smith 
77647c6ae99SBarry Smith    Notes:
777aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
778aa219208SBarry Smith    standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
77947c6ae99SBarry Smith    the standard 9-pt stencil.
78047c6ae99SBarry Smith 
781aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
782564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
783564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
78447c6ae99SBarry Smith 
785897f7067SBarry Smith    You must call DMSetUp() after this call before using this DM.
786897f7067SBarry Smith 
787897f7067SBarry Smith    If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call
788897f7067SBarry Smith    but before DMSetUp().
789897f7067SBarry Smith 
790aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
79199f0b487SRichard Tran Mills           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(),
7923bd220d7SPatrick Sanan           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges(),
7933bd220d7SPatrick Sanan           DMStagCreate2d()
79447c6ae99SBarry Smith 
79547c6ae99SBarry Smith @*/
796fe16a2e9SBarry Smith 
797bff4a2f0SMatthew G. Knepley PetscErrorCode  DMDACreate2d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMDAStencilType stencil_type,
7989a42bb27SBarry Smith                              PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
79947c6ae99SBarry Smith {
80047c6ae99SBarry Smith   PetscErrorCode ierr;
80147c6ae99SBarry Smith 
80247c6ae99SBarry Smith   PetscFunctionBegin;
803aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
804c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*da, 2);CHKERRQ(ierr);
805aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, 1);CHKERRQ(ierr);
806aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr);
807bff4a2f0SMatthew G. Knepley   ierr = DMDASetBoundaryType(*da, bx, by, DM_BOUNDARY_NONE);CHKERRQ(ierr);
808aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
809aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
810aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
8110298fd71SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, NULL);CHKERRQ(ierr);
81247c6ae99SBarry Smith   PetscFunctionReturn(0);
81347c6ae99SBarry Smith }
814