xref: /petsc/src/dm/impls/da/da2.c (revision 9566063d113dddea24716c546802770db7481bc0)
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   PetscMPIInt    rank;
8c9493c35SStefano Zampini   PetscBool      iascii,isdraw,isglvis,isbinary;
947c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
109a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
119a42bb27SBarry Smith   PetscBool ismatlab;
129a42bb27SBarry Smith #endif
1347c6ae99SBarry Smith 
1447c6ae99SBarry Smith   PetscFunctionBegin;
15*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank));
1647c6ae99SBarry Smith 
17*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
18*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw));
19*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis));
20*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary));
219a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
22*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab));
239a42bb27SBarry Smith #endif
2447c6ae99SBarry Smith   if (iascii) {
2547c6ae99SBarry Smith     PetscViewerFormat format;
2647c6ae99SBarry Smith 
27*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
2876a8abe0SBarry Smith     if (format == PETSC_VIEWER_LOAD_BALANCE) {
2976a8abe0SBarry Smith       PetscInt      i,nmax = 0,nmin = PETSC_MAX_INT,navg = 0,*nz,nzlocal;
3076a8abe0SBarry Smith       DMDALocalInfo info;
3176a8abe0SBarry Smith       PetscMPIInt   size;
32*9566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)da),&size));
33*9566063dSJacob Faibussowitsch       PetscCall(DMDAGetLocalInfo(da,&info));
3476a8abe0SBarry Smith       nzlocal = info.xm*info.ym;
35*9566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(size,&nz));
36*9566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Allgather(&nzlocal,1,MPIU_INT,nz,1,MPIU_INT,PetscObjectComm((PetscObject)da)));
3776a8abe0SBarry Smith       for (i=0; i<(PetscInt)size; i++) {
3876a8abe0SBarry Smith         nmax = PetscMax(nmax,nz[i]);
3976a8abe0SBarry Smith         nmin = PetscMin(nmin,nz[i]);
4076a8abe0SBarry Smith         navg += nz[i];
4176a8abe0SBarry Smith       }
42*9566063dSJacob Faibussowitsch       PetscCall(PetscFree(nz));
4376a8abe0SBarry Smith       navg = navg/size;
44*9566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  Load Balance - Grid Points: Min %D  avg %D  max %D\n",nmin,navg,nmax));
4576a8abe0SBarry Smith       PetscFunctionReturn(0);
4676a8abe0SBarry Smith     }
478ec8862eSJed Brown     if (format != PETSC_VIEWER_ASCII_VTK_DEPRECATED && format != PETSC_VIEWER_ASCII_VTK_CELL_DEPRECATED && format != PETSC_VIEWER_ASCII_GLVIS) {
48aa219208SBarry Smith       DMDALocalInfo info;
49*9566063dSJacob Faibussowitsch       PetscCall(DMDAGetLocalInfo(da,&info));
50*9566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
51*9566063dSJacob Faibussowitsch       PetscCall(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));
52*9566063dSJacob Faibussowitsch       PetscCall(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));
53*9566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
54*9566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
558135c375SStefano Zampini     } else if (format == PETSC_VIEWER_ASCII_GLVIS) {
56*9566063dSJacob Faibussowitsch       PetscCall(DMView_DA_GLVis(da,viewer));
573da9ae13SJed Brown     } else {
58*9566063dSJacob Faibussowitsch       PetscCall(DMView_DA_VTK(da,viewer));
5947c6ae99SBarry Smith     }
6047c6ae99SBarry Smith   } else if (isdraw) {
6147c6ae99SBarry Smith     PetscDraw      draw;
6247c6ae99SBarry Smith     double         ymin = -1*dd->s-1,ymax = dd->N+dd->s;
6347c6ae99SBarry Smith     double         xmin = -1*dd->s-1,xmax = dd->M+dd->s;
6447c6ae99SBarry Smith     double         x,y;
658ea3bf28SBarry Smith     PetscInt       base;
668ea3bf28SBarry Smith     const PetscInt *idx;
6747c6ae99SBarry Smith     char           node[10];
6847c6ae99SBarry Smith     PetscBool      isnull;
695f80ce2aSJacob Faibussowitsch     PetscErrorCode ierr;
7047c6ae99SBarry Smith 
71*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw));
72*9566063dSJacob Faibussowitsch     PetscCall(PetscDrawIsNull(draw,&isnull));
735b399a63SLisandro Dalcin     if (isnull) PetscFunctionReturn(0);
7447c6ae99SBarry Smith 
75*9566063dSJacob Faibussowitsch     PetscCall(PetscDrawCheckResizedWindow(draw));
76*9566063dSJacob Faibussowitsch     PetscCall(PetscDrawClear(draw));
77*9566063dSJacob Faibussowitsch     PetscCall(PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax));
785b399a63SLisandro Dalcin 
79*9566063dSJacob Faibussowitsch     ierr = PetscDrawCollectiveBegin(draw);PetscCall(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++) {
84*9566063dSJacob Faibussowitsch         PetscCall(PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK));
8547c6ae99SBarry Smith       }
8647c6ae99SBarry Smith       xmin = 0.0; xmax = dd->M - 1;
8747c6ae99SBarry Smith       for (ymin=0; ymin<dd->N; ymin++) {
88*9566063dSJacob Faibussowitsch         PetscCall(PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK));
8947c6ae99SBarry Smith       }
9047c6ae99SBarry Smith     }
91*9566063dSJacob Faibussowitsch     ierr = PetscDrawCollectiveEnd(draw);PetscCall(ierr);
92*9566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
93*9566063dSJacob Faibussowitsch     PetscCall(PetscDrawPause(draw));
9447c6ae99SBarry Smith 
95*9566063dSJacob Faibussowitsch     ierr = PetscDrawCollectiveBegin(draw);PetscCall(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;
98*9566063dSJacob Faibussowitsch     PetscCall(PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED));
99*9566063dSJacob Faibussowitsch     PetscCall(PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED));
100*9566063dSJacob Faibussowitsch     PetscCall(PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED));
101*9566063dSJacob Faibussowitsch     PetscCall(PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED));
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++) {
106*9566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(node,sizeof(node),"%d",(int)base++));
107*9566063dSJacob Faibussowitsch         PetscCall(PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node));
10847c6ae99SBarry Smith       }
10947c6ae99SBarry Smith     }
110*9566063dSJacob Faibussowitsch     ierr = PetscDrawCollectiveEnd(draw);PetscCall(ierr);
111*9566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
112*9566063dSJacob Faibussowitsch     PetscCall(PetscDrawPause(draw));
11347c6ae99SBarry Smith 
114*9566063dSJacob Faibussowitsch     ierr = PetscDrawCollectiveBegin(draw);PetscCall(ierr);
1155b399a63SLisandro Dalcin     /* overlay ghost numbers, useful for error checking */
116*9566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx));
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) {
121*9566063dSJacob Faibussowitsch           PetscCall(PetscSNPrintf(node,sizeof(node),"%d",(int)(idx[base/dd->w])));
122*9566063dSJacob Faibussowitsch           PetscCall(PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node));
12347c6ae99SBarry Smith         }
12447c6ae99SBarry Smith         base++;
12547c6ae99SBarry Smith       }
12647c6ae99SBarry Smith     }
127*9566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx));
128*9566063dSJacob Faibussowitsch     ierr = PetscDrawCollectiveEnd(draw);PetscCall(ierr);
129*9566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
130*9566063dSJacob Faibussowitsch     PetscCall(PetscDrawPause(draw));
131*9566063dSJacob Faibussowitsch     PetscCall(PetscDrawSave(draw));
132c9493c35SStefano Zampini   } else if (isglvis) {
133*9566063dSJacob Faibussowitsch     PetscCall(DMView_DA_GLVis(da,viewer));
1349a42bb27SBarry Smith   } else if (isbinary) {
135*9566063dSJacob Faibussowitsch     PetscCall(DMView_DA_Binary(da,viewer));
1369a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1379a42bb27SBarry Smith   } else if (ismatlab) {
138*9566063dSJacob Faibussowitsch     PetscCall(DMView_DA_Matlab(da,viewer));
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   PetscInt       gI,nI;
15947c6ae99SBarry Smith   MatStencil     stencil;
160aa219208SBarry Smith   DMDALocalInfo  info;
16147c6ae99SBarry Smith 
16247c6ae99SBarry Smith   PetscFunctionBegin;
163*9566063dSJacob Faibussowitsch   PetscCall((*ctx->func)(0,U,a,ctx->funcctx));
164*9566063dSJacob Faibussowitsch   PetscCall((*ctx->funcisetbase)(U,ctx->funcctx));
16547c6ae99SBarry Smith 
166*9566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U,&ww));
167*9566063dSJacob Faibussowitsch   PetscCall(VecGetArray(a,&aa));
16847c6ae99SBarry Smith 
16947c6ae99SBarry Smith   nI = 0;
17047c6ae99SBarry Smith   h  = ww[gI];
17147c6ae99SBarry Smith   if (h == 0.0) h = 1.0;
17247c6ae99SBarry Smith   if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin;
17347c6ae99SBarry Smith   else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
17447c6ae99SBarry Smith   h *= epsilon;
17547c6ae99SBarry Smith 
17647c6ae99SBarry Smith   ww[gI] += h;
177*9566063dSJacob Faibussowitsch   PetscCall((*ctx->funci)(i,w,&v,ctx->funcctx));
17847c6ae99SBarry Smith   aa[nI]  = (v - aa[nI])/h;
17947c6ae99SBarry Smith   ww[gI] -= h;
18047c6ae99SBarry Smith   nI++;
1818865f1eaSKarl Rupp 
182*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U,&ww));
183*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(a,&aa));
18447c6ae99SBarry Smith   PetscFunctionReturn(0);
18547c6ae99SBarry Smith }
18647c6ae99SBarry Smith #endif
18747c6ae99SBarry Smith 
1887087cfbeSBarry Smith PetscErrorCode  DMSetUp_DA_2D(DM da)
18947c6ae99SBarry Smith {
19047c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
19147c6ae99SBarry Smith   const PetscInt   M            = dd->M;
19247c6ae99SBarry Smith   const PetscInt   N            = dd->N;
19347c6ae99SBarry Smith   PetscInt         m            = dd->m;
19447c6ae99SBarry Smith   PetscInt         n            = dd->n;
19547c6ae99SBarry Smith   const PetscInt   dof          = dd->w;
19647c6ae99SBarry Smith   const PetscInt   s            = dd->s;
197bff4a2f0SMatthew G. Knepley   DMBoundaryType   bx           = dd->bx;
198bff4a2f0SMatthew G. Knepley   DMBoundaryType   by           = dd->by;
19919fd82e9SBarry Smith   DMDAStencilType  stencil_type = dd->stencil_type;
20047c6ae99SBarry Smith   PetscInt         *lx          = dd->lx;
20147c6ae99SBarry Smith   PetscInt         *ly          = dd->ly;
20247c6ae99SBarry Smith   MPI_Comm         comm;
20347c6ae99SBarry Smith   PetscMPIInt      rank,size;
204bd1fc5aeSBarry Smith   PetscInt         xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,IXs,IXe,IYs,IYe;
2058ea3bf28SBarry Smith   PetscInt         up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn;
206db87c5ecSEthan Coon   PetscInt         xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
20747c6ae99SBarry Smith   PetscInt         s_x,s_y; /* s proportionalized to w */
20847c6ae99SBarry Smith   PetscInt         sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
20947c6ae99SBarry Smith   Vec              local,global;
210bd1fc5aeSBarry Smith   VecScatter       gtol;
21145b6f7e9SBarry Smith   IS               to,from;
21247c6ae99SBarry Smith 
21347c6ae99SBarry Smith   PetscFunctionBegin;
2147a8be351SBarry 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");
215*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)da,&comm));
2163855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES)
2177a8be351SBarry 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);
2183855c12bSBarry Smith #endif
21947c6ae99SBarry Smith 
220*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
221*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
22247c6ae99SBarry Smith 
2237d310018SBarry Smith   dd->p = 1;
22447c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
2257a8be351SBarry Smith     PetscCheck(m >= 1,comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
2267a8be351SBarry Smith     else PetscCheck(m <= size,comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
22747c6ae99SBarry Smith   }
22847c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
2297a8be351SBarry Smith     PetscCheck(n >= 1,comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
2307a8be351SBarry Smith     else PetscCheck(n <= size,comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
23147c6ae99SBarry Smith   }
23247c6ae99SBarry Smith 
23347c6ae99SBarry Smith   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
23447c6ae99SBarry Smith     if (n != PETSC_DECIDE) {
23547c6ae99SBarry Smith       m = size/n;
23647c6ae99SBarry Smith     } else if (m != PETSC_DECIDE) {
23747c6ae99SBarry Smith       n = size/m;
23847c6ae99SBarry Smith     } else {
23947c6ae99SBarry Smith       /* try for squarish distribution */
240369cc0aeSBarry Smith       m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N)));
24147c6ae99SBarry Smith       if (!m) m = 1;
24247c6ae99SBarry Smith       while (m > 0) {
24347c6ae99SBarry Smith         n = size/m;
24447c6ae99SBarry Smith         if (m*n == size) break;
24547c6ae99SBarry Smith         m--;
24647c6ae99SBarry Smith       }
24747c6ae99SBarry Smith       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
24847c6ae99SBarry Smith     }
2497a8be351SBarry Smith     PetscCheck(m*n == size,comm,PETSC_ERR_PLIB,"Unable to create partition, check the size of the communicator and input m and n ");
2507a8be351SBarry Smith   } else PetscCheck(m*n == size,comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
25147c6ae99SBarry Smith 
2527a8be351SBarry Smith   PetscCheck(M >= m,comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
2537a8be351SBarry Smith   PetscCheck(N >= n,comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
25447c6ae99SBarry Smith 
25547c6ae99SBarry Smith   /*
25647c6ae99SBarry Smith      Determine locally owned region
25747c6ae99SBarry Smith      xs is the first local node number, x is the number of local nodes
25847c6ae99SBarry Smith   */
25947c6ae99SBarry Smith   if (!lx) {
260*9566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m, &dd->lx));
26147c6ae99SBarry Smith     lx   = dd->lx;
26247c6ae99SBarry Smith     for (i=0; i<m; i++) {
26347c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > i);
26447c6ae99SBarry Smith     }
26547c6ae99SBarry Smith   }
26647c6ae99SBarry Smith   x  = lx[rank % m];
26747c6ae99SBarry Smith   xs = 0;
26847c6ae99SBarry Smith   for (i=0; i<(rank % m); i++) {
26947c6ae99SBarry Smith     xs += lx[i];
27047c6ae99SBarry Smith   }
27176bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
27247c6ae99SBarry Smith     left = xs;
27347c6ae99SBarry Smith     for (i=(rank % m); i<m; i++) {
27447c6ae99SBarry Smith       left += lx[i];
27547c6ae99SBarry Smith     }
2767a8be351SBarry Smith     PetscCheck(left == M,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M: %D %D",left,M);
27776bd3646SJed Brown   }
27847c6ae99SBarry Smith 
27947c6ae99SBarry Smith   /*
28047c6ae99SBarry Smith      Determine locally owned region
28147c6ae99SBarry Smith      ys is the first local node number, y is the number of local nodes
28247c6ae99SBarry Smith   */
28347c6ae99SBarry Smith   if (!ly) {
284*9566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &dd->ly));
28547c6ae99SBarry Smith     ly   = dd->ly;
28647c6ae99SBarry Smith     for (i=0; i<n; i++) {
28747c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > i);
28847c6ae99SBarry Smith     }
28947c6ae99SBarry Smith   }
29047c6ae99SBarry Smith   y  = ly[rank/m];
29147c6ae99SBarry Smith   ys = 0;
29247c6ae99SBarry Smith   for (i=0; i<(rank/m); i++) {
29347c6ae99SBarry Smith     ys += ly[i];
29447c6ae99SBarry Smith   }
29576bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
29647c6ae99SBarry Smith     left = ys;
29747c6ae99SBarry Smith     for (i=(rank/m); i<n; i++) {
29847c6ae99SBarry Smith       left += ly[i];
29947c6ae99SBarry Smith     }
3007a8be351SBarry Smith     PetscCheck(left == N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of ly across processors not equal to N: %D %D",left,N);
30176bd3646SJed Brown   }
30247c6ae99SBarry Smith 
303bcea557cSEthan Coon   /*
304bcea557cSEthan Coon    check if the scatter requires more than one process neighbor or wraps around
305bcea557cSEthan Coon    the domain more than once
306bcea557cSEthan Coon   */
3077a8be351SBarry 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);
3087a8be351SBarry 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);
30947c6ae99SBarry Smith   xe = xs + x;
31047c6ae99SBarry Smith   ye = ys + y;
31147c6ae99SBarry Smith 
312ce00eea3SSatish Balay   /* determine ghost region (Xs) and region scattered into (IXs)  */
313d9c9ebe5SPeter Brune   if (xs-s > 0) {
314d9c9ebe5SPeter Brune     Xs = xs - s; IXs = xs - s;
31588661749SPeter Brune   } else {
31688661749SPeter Brune     if (bx) {
31788661749SPeter Brune       Xs = xs - s;
31888661749SPeter Brune     } else {
31988661749SPeter Brune       Xs = 0;
32088661749SPeter Brune     }
32188661749SPeter Brune     IXs = 0;
32288661749SPeter Brune   }
323d9c9ebe5SPeter Brune   if (xe+s <= M) {
324d9c9ebe5SPeter Brune     Xe = xe + s; IXe = xe + s;
32588661749SPeter Brune   } else {
32688661749SPeter Brune     if (bx) {
327d9c9ebe5SPeter Brune       Xs = xs - s; Xe = xe + s;
32888661749SPeter Brune     } else {
32988661749SPeter Brune       Xe = M;
33088661749SPeter Brune     }
33188661749SPeter Brune     IXe = M;
33288661749SPeter Brune   }
33347c6ae99SBarry Smith 
334bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) {
335d9c9ebe5SPeter Brune     IXs = xs - s;
336d9c9ebe5SPeter Brune     IXe = xe + s;
337d9c9ebe5SPeter Brune     Xs  = xs - s;
338d9c9ebe5SPeter Brune     Xe  = xe + s;
33988661749SPeter Brune   }
34047c6ae99SBarry Smith 
341d9c9ebe5SPeter Brune   if (ys-s > 0) {
342d9c9ebe5SPeter Brune     Ys = ys - s; IYs = ys - s;
34388661749SPeter Brune   } else {
34488661749SPeter Brune     if (by) {
34588661749SPeter Brune       Ys = ys - s;
34688661749SPeter Brune     } else {
34788661749SPeter Brune       Ys = 0;
34888661749SPeter Brune     }
34988661749SPeter Brune     IYs = 0;
35088661749SPeter Brune   }
351d9c9ebe5SPeter Brune   if (ye+s <= N) {
352d9c9ebe5SPeter Brune     Ye = ye + s; IYe = ye + s;
35388661749SPeter Brune   } else {
35488661749SPeter Brune     if (by) {
35588661749SPeter Brune       Ye = ye + s;
35688661749SPeter Brune     } else {
35788661749SPeter Brune       Ye = N;
35888661749SPeter Brune     }
35988661749SPeter Brune     IYe = N;
36088661749SPeter Brune   }
36188661749SPeter Brune 
362bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) {
363d9c9ebe5SPeter Brune     IYs = ys - s;
364d9c9ebe5SPeter Brune     IYe = ye + s;
365d9c9ebe5SPeter Brune     Ys  = ys - s;
366d9c9ebe5SPeter Brune     Ye  = ye + s;
36788661749SPeter Brune   }
36888661749SPeter Brune 
36988661749SPeter Brune   /* stencil length in each direction */
370d9c9ebe5SPeter Brune   s_x = s;
371d9c9ebe5SPeter Brune   s_y = s;
37247c6ae99SBarry Smith 
37347c6ae99SBarry Smith   /* determine starting point of each processor */
37447c6ae99SBarry Smith   nn       = x*y;
375*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(size+1,&bases,size,&ldims));
376*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm));
37747c6ae99SBarry Smith   bases[0] = 0;
37847c6ae99SBarry Smith   for (i=1; i<=size; i++) {
37947c6ae99SBarry Smith     bases[i] = ldims[i-1];
38047c6ae99SBarry Smith   }
38147c6ae99SBarry Smith   for (i=1; i<=size; i++) {
38247c6ae99SBarry Smith     bases[i] += bases[i-1];
38347c6ae99SBarry Smith   }
384ce00eea3SSatish Balay   base = bases[rank]*dof;
38547c6ae99SBarry Smith 
38647c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
387ce00eea3SSatish Balay   dd->Nlocal = x*y*dof;
388*9566063dSJacob Faibussowitsch   PetscCall(VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global));
389ce00eea3SSatish Balay   dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof;
390*9566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local));
39147c6ae99SBarry Smith 
3924104a7a0SPatrick Sanan   /* generate global to local vector scatter and local to global mapping*/
39347c6ae99SBarry Smith 
394ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
395ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
396*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((IXe-IXs)*(IYe-IYs),&idx));
397d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_BOX) {
398ce00eea3SSatish Balay     left  = IXs - Xs; right = left + (IXe-IXs);
399ce00eea3SSatish Balay     down  = IYs - Ys; up = down + (IYe-IYs);
400ce00eea3SSatish Balay     count = 0;
401ce00eea3SSatish Balay     for (i=down; i<up; i++) {
402ce00eea3SSatish Balay       for (j=left; j<right; j++) {
403ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
404ce00eea3SSatish Balay       }
405ce00eea3SSatish Balay     }
406*9566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to));
407ce00eea3SSatish Balay 
40847c6ae99SBarry Smith   } else {
40947c6ae99SBarry Smith     /* must drop into cross shape region */
41047c6ae99SBarry Smith     /*       ---------|
41147c6ae99SBarry Smith             |  top    |
412ce00eea3SSatish Balay          |---         ---| up
41347c6ae99SBarry Smith          |   middle      |
41447c6ae99SBarry Smith          |               |
415ce00eea3SSatish Balay          ----         ---- down
41647c6ae99SBarry Smith             | bottom  |
41747c6ae99SBarry Smith             -----------
41847c6ae99SBarry Smith          Xs xs        xe Xe */
419ce00eea3SSatish Balay     left  = xs - Xs; right = left + x;
420ce00eea3SSatish Balay     down  = ys - Ys; up = down + y;
42147c6ae99SBarry Smith     count = 0;
422ce00eea3SSatish Balay     /* bottom */
423ce00eea3SSatish Balay     for (i=(IYs-Ys); i<down; i++) {
424ce00eea3SSatish Balay       for (j=left; j<right; j++) {
425ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
42647c6ae99SBarry Smith       }
42747c6ae99SBarry Smith     }
42847c6ae99SBarry Smith     /* middle */
42947c6ae99SBarry Smith     for (i=down; i<up; i++) {
430ce00eea3SSatish Balay       for (j=(IXs-Xs); j<(IXe-Xs); j++) {
431ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
43247c6ae99SBarry Smith       }
43347c6ae99SBarry Smith     }
43447c6ae99SBarry Smith     /* top */
435ce00eea3SSatish Balay     for (i=up; i<up+IYe-ye; i++) {
436ce00eea3SSatish Balay       for (j=left; j<right; j++) {
437ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
43847c6ae99SBarry Smith       }
43947c6ae99SBarry Smith     }
440*9566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to));
44147c6ae99SBarry Smith   }
44247c6ae99SBarry Smith 
44347c6ae99SBarry Smith   /* determine who lies on each side of us stored in    n6 n7 n8
44447c6ae99SBarry Smith                                                         n3    n5
44547c6ae99SBarry Smith                                                         n0 n1 n2
44647c6ae99SBarry Smith   */
44747c6ae99SBarry Smith 
44847c6ae99SBarry Smith   /* Assume the Non-Periodic Case */
44947c6ae99SBarry Smith   n1 = rank - m;
45047c6ae99SBarry Smith   if (rank % m) {
45147c6ae99SBarry Smith     n0 = n1 - 1;
45247c6ae99SBarry Smith   } else {
45347c6ae99SBarry Smith     n0 = -1;
45447c6ae99SBarry Smith   }
45547c6ae99SBarry Smith   if ((rank+1) % m) {
45647c6ae99SBarry Smith     n2 = n1 + 1;
45747c6ae99SBarry Smith     n5 = rank + 1;
45847c6ae99SBarry Smith     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
45947c6ae99SBarry Smith   } else {
46047c6ae99SBarry Smith     n2 = -1; n5 = -1; n8 = -1;
46147c6ae99SBarry Smith   }
46247c6ae99SBarry Smith   if (rank % m) {
46347c6ae99SBarry Smith     n3 = rank - 1;
46447c6ae99SBarry Smith     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
46547c6ae99SBarry Smith   } else {
46647c6ae99SBarry Smith     n3 = -1; n6 = -1;
46747c6ae99SBarry Smith   }
46847c6ae99SBarry Smith   n7 = rank + m; if (n7 >= m*n) n7 = -1;
46947c6ae99SBarry Smith 
470bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC && by == DM_BOUNDARY_PERIODIC) {
47147c6ae99SBarry Smith     /* Modify for Periodic Cases */
47247c6ae99SBarry Smith     /* Handle all four corners */
47347c6ae99SBarry Smith     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
47447c6ae99SBarry Smith     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
47547c6ae99SBarry Smith     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
47647c6ae99SBarry Smith     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
47747c6ae99SBarry Smith 
47847c6ae99SBarry Smith     /* Handle Top and Bottom Sides */
47947c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
48047c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
48147c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
48247c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
48347c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
48447c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
48547c6ae99SBarry Smith 
48647c6ae99SBarry Smith     /* Handle Left and Right Sides */
48747c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
48847c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
48947c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
49047c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
49147c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
49247c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
493bff4a2f0SMatthew G. Knepley   } else if (by == DM_BOUNDARY_PERIODIC) {  /* Handle Top and Bottom Sides */
494ce00eea3SSatish Balay     if (n1 < 0) n1 = rank + m * (n-1);
495ce00eea3SSatish Balay     if (n7 < 0) n7 = rank - m * (n-1);
496ce00eea3SSatish Balay     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
497ce00eea3SSatish Balay     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
498ce00eea3SSatish Balay     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
499ce00eea3SSatish Balay     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
500bff4a2f0SMatthew G. Knepley   } else if (bx == DM_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */
501ce00eea3SSatish Balay     if (n3 < 0) n3 = rank + (m-1);
502ce00eea3SSatish Balay     if (n5 < 0) n5 = rank - (m-1);
503ce00eea3SSatish Balay     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
504ce00eea3SSatish Balay     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
505ce00eea3SSatish Balay     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
506ce00eea3SSatish Balay     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
50747c6ae99SBarry Smith   }
508ce00eea3SSatish Balay 
509*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(9,&dd->neighbors));
5108865f1eaSKarl Rupp 
51147c6ae99SBarry Smith   dd->neighbors[0] = n0;
51247c6ae99SBarry Smith   dd->neighbors[1] = n1;
51347c6ae99SBarry Smith   dd->neighbors[2] = n2;
51447c6ae99SBarry Smith   dd->neighbors[3] = n3;
51547c6ae99SBarry Smith   dd->neighbors[4] = rank;
51647c6ae99SBarry Smith   dd->neighbors[5] = n5;
51747c6ae99SBarry Smith   dd->neighbors[6] = n6;
51847c6ae99SBarry Smith   dd->neighbors[7] = n7;
51947c6ae99SBarry Smith   dd->neighbors[8] = n8;
52047c6ae99SBarry Smith 
521d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
52247c6ae99SBarry Smith     /* save corner processor numbers */
52347c6ae99SBarry Smith     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
52447c6ae99SBarry Smith     n0  = n2 = n6 = n8 = -1;
52547c6ae99SBarry Smith   }
52647c6ae99SBarry Smith 
527*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((Xe-Xs)*(Ye-Ys),&idx));
52847c6ae99SBarry Smith 
529ce00eea3SSatish Balay   nn = 0;
53047c6ae99SBarry Smith   xbase = bases[rank];
53147c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
53247c6ae99SBarry Smith     if (n0 >= 0) { /* left below */
533ce00eea3SSatish Balay       x_t = lx[n0 % m];
53447c6ae99SBarry Smith       y_t = ly[(n0/m)];
53547c6ae99SBarry Smith       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
5368865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
53747c6ae99SBarry Smith     }
538ac119b13SBarry Smith 
53947c6ae99SBarry Smith     if (n1 >= 0) { /* directly below */
54047c6ae99SBarry Smith       x_t = x;
54147c6ae99SBarry Smith       y_t = ly[(n1/m)];
54247c6ae99SBarry Smith       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
5438865f1eaSKarl Rupp       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
544bff4a2f0SMatthew G. Knepley     } else if (by == DM_BOUNDARY_MIRROR) {
5458865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
54647c6ae99SBarry Smith     }
547ac119b13SBarry Smith 
54847c6ae99SBarry Smith     if (n2 >= 0) { /* right below */
549ce00eea3SSatish Balay       x_t = lx[n2 % m];
55047c6ae99SBarry Smith       y_t = ly[(n2/m)];
55147c6ae99SBarry Smith       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
5528865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
55347c6ae99SBarry Smith     }
55447c6ae99SBarry Smith   }
55547c6ae99SBarry Smith 
55647c6ae99SBarry Smith   for (i=0; i<y; i++) {
55747c6ae99SBarry Smith     if (n3 >= 0) { /* directly left */
558ce00eea3SSatish Balay       x_t = lx[n3 % m];
55947c6ae99SBarry Smith       /* y_t = y; */
56047c6ae99SBarry Smith       s_t = bases[n3] + (i+1)*x_t - s_x;
5618865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
562bff4a2f0SMatthew G. Knepley     } else if (bx == DM_BOUNDARY_MIRROR) {
5638865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
56447c6ae99SBarry Smith     }
56547c6ae99SBarry Smith 
5668865f1eaSKarl Rupp     for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
56747c6ae99SBarry Smith 
56847c6ae99SBarry Smith     if (n5 >= 0) { /* directly right */
569ce00eea3SSatish Balay       x_t = lx[n5 % m];
57047c6ae99SBarry Smith       /* y_t = y; */
57147c6ae99SBarry Smith       s_t = bases[n5] + (i)*x_t;
5728865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
573bff4a2f0SMatthew G. Knepley     } else if (bx == DM_BOUNDARY_MIRROR) {
5748865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
57547c6ae99SBarry Smith     }
57647c6ae99SBarry Smith   }
57747c6ae99SBarry Smith 
57847c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
57947c6ae99SBarry Smith     if (n6 >= 0) { /* left above */
580ce00eea3SSatish Balay       x_t = lx[n6 % m];
58147c6ae99SBarry Smith       /* y_t = ly[(n6/m)]; */
58247c6ae99SBarry Smith       s_t = bases[n6] + (i)*x_t - s_x;
5838865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
58447c6ae99SBarry Smith     }
585ac119b13SBarry Smith 
58647c6ae99SBarry Smith     if (n7 >= 0) { /* directly above */
58747c6ae99SBarry Smith       x_t = x;
58847c6ae99SBarry Smith       /* y_t = ly[(n7/m)]; */
58947c6ae99SBarry Smith       s_t = bases[n7] + (i-1)*x_t;
5908865f1eaSKarl Rupp       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
591bff4a2f0SMatthew G. Knepley     } else if (by == DM_BOUNDARY_MIRROR) {
5928865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
59347c6ae99SBarry Smith     }
594ac119b13SBarry Smith 
59547c6ae99SBarry Smith     if (n8 >= 0) { /* right above */
596ce00eea3SSatish Balay       x_t = lx[n8 % m];
59747c6ae99SBarry Smith       /* y_t = ly[(n8/m)]; */
59847c6ae99SBarry Smith       s_t = bases[n8] + (i-1)*x_t;
5998865f1eaSKarl Rupp       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
60047c6ae99SBarry Smith     }
60147c6ae99SBarry Smith   }
60247c6ae99SBarry Smith 
603*9566063dSJacob Faibussowitsch   PetscCall(ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from));
604*9566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(global,from,local,to,&gtol));
605*9566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)da,(PetscObject)gtol));
606*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
607*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
60847c6ae99SBarry Smith 
609d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
610ce00eea3SSatish Balay     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
611ce00eea3SSatish Balay   }
612ce00eea3SSatish Balay 
613288e7d53SBarry Smith   if (((stencil_type == DMDA_STENCIL_STAR)  || (bx && bx != DM_BOUNDARY_PERIODIC) || (by && by != DM_BOUNDARY_PERIODIC))) {
61447c6ae99SBarry Smith     /*
61547c6ae99SBarry Smith         Recompute the local to global mappings, this time keeping the
616ce00eea3SSatish Balay       information about the cross corner processor numbers and any ghosted
617ce00eea3SSatish Balay       but not periodic indices.
61847c6ae99SBarry Smith     */
61947c6ae99SBarry Smith     nn    = 0;
62047c6ae99SBarry Smith     xbase = bases[rank];
62147c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
62247c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
623ce00eea3SSatish Balay         x_t = lx[n0 % m];
62447c6ae99SBarry Smith         y_t = ly[(n0/m)];
62547c6ae99SBarry Smith         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
6268865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
627ce00eea3SSatish Balay       } else if (xs-Xs > 0 && ys-Ys > 0) {
6288865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
62947c6ae99SBarry Smith       }
63047c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
63147c6ae99SBarry Smith         x_t = x;
63247c6ae99SBarry Smith         y_t = ly[(n1/m)];
63347c6ae99SBarry Smith         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
6348865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
635ce00eea3SSatish Balay       } else if (ys-Ys > 0) {
636bff4a2f0SMatthew G. Knepley         if (by == DM_BOUNDARY_MIRROR) {
6378865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
638624904c4SBarry Smith         } else {
6398865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
64047c6ae99SBarry Smith         }
641624904c4SBarry Smith       }
64247c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
643ce00eea3SSatish Balay         x_t = lx[n2 % m];
64447c6ae99SBarry Smith         y_t = ly[(n2/m)];
64547c6ae99SBarry Smith         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
6468865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
647ce00eea3SSatish Balay       } else if (Xe-xe> 0 && ys-Ys > 0) {
6488865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
64947c6ae99SBarry Smith       }
65047c6ae99SBarry Smith     }
65147c6ae99SBarry Smith 
65247c6ae99SBarry Smith     for (i=0; i<y; i++) {
65347c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
654ce00eea3SSatish Balay         x_t = lx[n3 % m];
65547c6ae99SBarry Smith         /* y_t = y; */
65647c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x;
6578865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
658ce00eea3SSatish Balay       } else if (xs-Xs > 0) {
659bff4a2f0SMatthew G. Knepley         if (bx == DM_BOUNDARY_MIRROR) {
6608865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
661624904c4SBarry Smith         } else {
6628865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
66347c6ae99SBarry Smith         }
664624904c4SBarry Smith       }
66547c6ae99SBarry Smith 
6668865f1eaSKarl Rupp       for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
66747c6ae99SBarry Smith 
66847c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
669ce00eea3SSatish Balay         x_t = lx[n5 % m];
67047c6ae99SBarry Smith         /* y_t = y; */
67147c6ae99SBarry Smith         s_t = bases[n5] + (i)*x_t;
6728865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
673ce00eea3SSatish Balay       } else if (Xe-xe > 0) {
674bff4a2f0SMatthew G. Knepley         if (bx == DM_BOUNDARY_MIRROR) {
6758865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
676624904c4SBarry Smith         } else {
6778865f1eaSKarl Rupp           for (j=0; j<s_x; j++) idx[nn++] = -1;
67847c6ae99SBarry Smith         }
67947c6ae99SBarry Smith       }
680624904c4SBarry Smith     }
68147c6ae99SBarry Smith 
68247c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
68347c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
684ce00eea3SSatish Balay         x_t = lx[n6 % m];
68547c6ae99SBarry Smith         /* y_t = ly[(n6/m)]; */
68647c6ae99SBarry Smith         s_t = bases[n6] + (i)*x_t - s_x;
6878865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
688ce00eea3SSatish Balay       } else if (xs-Xs > 0 && Ye-ye > 0) {
6898865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
69047c6ae99SBarry Smith       }
69147c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
69247c6ae99SBarry Smith         x_t = x;
69347c6ae99SBarry Smith         /* y_t = ly[(n7/m)]; */
69447c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t;
6958865f1eaSKarl Rupp         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
696ce00eea3SSatish Balay       } else if (Ye-ye > 0) {
697bff4a2f0SMatthew G. Knepley         if (by == DM_BOUNDARY_MIRROR) {
6988865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
699624904c4SBarry Smith         } else {
7008865f1eaSKarl Rupp           for (j=0; j<x; j++) idx[nn++] = -1;
70147c6ae99SBarry Smith         }
702624904c4SBarry Smith       }
70347c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
704ce00eea3SSatish Balay         x_t = lx[n8 % m];
70547c6ae99SBarry Smith         /* y_t = ly[(n8/m)]; */
70647c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t;
7078865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
708ce00eea3SSatish Balay       } else if (Xe-xe > 0 && Ye-ye > 0) {
7098865f1eaSKarl Rupp         for (j=0; j<s_x; j++) idx[nn++] = -1;
71047c6ae99SBarry Smith       }
71147c6ae99SBarry Smith     }
71247c6ae99SBarry Smith   }
713ce00eea3SSatish Balay   /*
714ce00eea3SSatish Balay      Set the local to global ordering in the global vector, this allows use
715ce00eea3SSatish Balay      of VecSetValuesLocal().
716ce00eea3SSatish Balay   */
717*9566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap));
718*9566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap));
71947c6ae99SBarry Smith 
720*9566063dSJacob Faibussowitsch   PetscCall(PetscFree2(bases,ldims));
72147c6ae99SBarry Smith   dd->m = m;  dd->n  = n;
722ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
723ce00eea3SSatish Balay   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
724ce00eea3SSatish Balay   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
72547c6ae99SBarry Smith 
726*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&local));
727*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global));
72847c6ae99SBarry Smith 
72947c6ae99SBarry Smith   dd->gtol      = gtol;
73047c6ae99SBarry Smith   dd->base      = base;
7319a42bb27SBarry Smith   da->ops->view = DMView_DA_2d;
7320298fd71SBarry Smith   dd->ltol      = NULL;
7330298fd71SBarry Smith   dd->ao        = NULL;
73447c6ae99SBarry Smith   PetscFunctionReturn(0);
73547c6ae99SBarry Smith }
73647c6ae99SBarry Smith 
73747c6ae99SBarry Smith /*@C
738aa219208SBarry Smith    DMDACreate2d -  Creates an object that will manage the communication of  two-dimensional
73947c6ae99SBarry Smith    regular array data that is distributed across some processors.
74047c6ae99SBarry Smith 
741d083f849SBarry Smith    Collective
74247c6ae99SBarry Smith 
74347c6ae99SBarry Smith    Input Parameters:
74447c6ae99SBarry Smith +  comm - MPI communicator
7451321219cSEthan Coon .  bx,by - type of ghost nodes the array have.
746bff4a2f0SMatthew G. Knepley          Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC.
747aa219208SBarry Smith .  stencil_type - stencil type.  Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
748897f7067SBarry Smith .  M,N - global dimension in each direction of the array
74947c6ae99SBarry Smith .  m,n - corresponding number of processors in each dimension
75047c6ae99SBarry Smith          (or PETSC_DECIDE to have calculated)
75147c6ae99SBarry Smith .  dof - number of degrees of freedom per node
75247c6ae99SBarry Smith .  s - stencil width
75347c6ae99SBarry Smith -  lx, ly - arrays containing the number of nodes in each cell along
7540298fd71SBarry Smith            the x and y coordinates, or NULL. If non-null, these
75547c6ae99SBarry Smith            must be of length as m and n, and the corresponding
75647c6ae99SBarry Smith            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
75747c6ae99SBarry Smith            must be M, and the sum of the ly[] entries must be N.
75847c6ae99SBarry Smith 
75947c6ae99SBarry Smith    Output Parameter:
76047c6ae99SBarry Smith .  da - the resulting distributed array object
76147c6ae99SBarry Smith 
76247c6ae99SBarry Smith    Options Database Key:
763706a11cbSBarry Smith +  -dm_view - Calls DMView() at the conclusion of DMDACreate2d()
764e43dc8daSBarry Smith .  -da_grid_x <nx> - number of grid points in x direction
765e43dc8daSBarry Smith .  -da_grid_y <ny> - number of grid points in y direction
76647c6ae99SBarry Smith .  -da_processors_x <nx> - number of processors in x direction
76747c6ae99SBarry Smith .  -da_processors_y <ny> - number of processors in y direction
768e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement ratio in x direction
769e0f5d30fSBarry Smith .  -da_refine_y <ry> - refinement ratio in y direction
770e43dc8daSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating
771e0f5d30fSBarry Smith 
77247c6ae99SBarry Smith    Level: beginner
77347c6ae99SBarry Smith 
77447c6ae99SBarry Smith    Notes:
775aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
776aa219208SBarry Smith    standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
77747c6ae99SBarry Smith    the standard 9-pt stencil.
77847c6ae99SBarry Smith 
779aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
780564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
781564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
78247c6ae99SBarry Smith 
783897f7067SBarry Smith    You must call DMSetUp() after this call before using this DM.
784897f7067SBarry Smith 
785897f7067SBarry Smith    If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call
786897f7067SBarry Smith    but before DMSetUp().
787897f7067SBarry Smith 
788aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
78999f0b487SRichard Tran Mills           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(),
7903bd220d7SPatrick Sanan           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges(),
7913bd220d7SPatrick Sanan           DMStagCreate2d()
79247c6ae99SBarry Smith 
79347c6ae99SBarry Smith @*/
794fe16a2e9SBarry Smith 
795bff4a2f0SMatthew G. Knepley PetscErrorCode  DMDACreate2d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMDAStencilType stencil_type,
7969a42bb27SBarry Smith                              PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
79747c6ae99SBarry Smith {
79847c6ae99SBarry Smith   PetscFunctionBegin;
799*9566063dSJacob Faibussowitsch   PetscCall(DMDACreate(comm, da));
800*9566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*da, 2));
801*9566063dSJacob Faibussowitsch   PetscCall(DMDASetSizes(*da, M, N, 1));
802*9566063dSJacob Faibussowitsch   PetscCall(DMDASetNumProcs(*da, m, n, PETSC_DECIDE));
803*9566063dSJacob Faibussowitsch   PetscCall(DMDASetBoundaryType(*da, bx, by, DM_BOUNDARY_NONE));
804*9566063dSJacob Faibussowitsch   PetscCall(DMDASetDof(*da, dof));
805*9566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilType(*da, stencil_type));
806*9566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilWidth(*da, s));
807*9566063dSJacob Faibussowitsch   PetscCall(DMDASetOwnershipRanges(*da, lx, ly, NULL));
80847c6ae99SBarry Smith   PetscFunctionReturn(0);
80947c6ae99SBarry Smith }
810