xref: /petsc/src/dm/impls/da/da2.c (revision c6228bba3efb3955e5e29a4ef79e8b65616d20d2)
1 
2 #include <petsc/private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
3 #include <petscdraw.h>
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "DMView_DA_2d"
7 PetscErrorCode DMView_DA_2d(DM da,PetscViewer viewer)
8 {
9   PetscErrorCode ierr;
10   PetscMPIInt    rank;
11   PetscBool      iascii,isdraw,isbinary;
12   DM_DA          *dd = (DM_DA*)da->data;
13 #if defined(PETSC_HAVE_MATLAB_ENGINE)
14   PetscBool ismatlab;
15 #endif
16 
17   PetscFunctionBegin;
18   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr);
19 
20   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
21   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
22   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
23 #if defined(PETSC_HAVE_MATLAB_ENGINE)
24   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr);
25 #endif
26   if (iascii) {
27     PetscViewerFormat format;
28 
29     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
30     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) {
31       DMDALocalInfo info;
32       ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
33       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
34       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);
35       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);
36       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
37       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
38     } else {
39       ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr);
40     }
41   } else if (isdraw) {
42     PetscDraw      draw;
43     double         ymin = -1*dd->s-1,ymax = dd->N+dd->s;
44     double         xmin = -1*dd->s-1,xmax = dd->M+dd->s;
45     double         x,y;
46     PetscInt       base;
47     const PetscInt *idx;
48     char           node[10];
49     PetscBool      isnull;
50 
51     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
52     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
53     if (isnull) PetscFunctionReturn(0);
54 
55     ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
56     ierr = PetscDrawClear(draw);CHKERRQ(ierr);
57 
58     ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
59 
60     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
61     /* first processor draw all node lines */
62     if (!rank) {
63       ymin = 0.0; ymax = dd->N - 1;
64       for (xmin=0; xmin<dd->M; xmin++) {
65         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
66       }
67       xmin = 0.0; xmax = dd->M - 1;
68       for (ymin=0; ymin<dd->N; ymin++) {
69         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
70       }
71     }
72     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
73 
74     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
75     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
76 
77     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
78     /* draw my box */
79     xmin = dd->xs/dd->w; xmax =(dd->xe-1)/dd->w; ymin = dd->ys; ymax = dd->ye - 1;
80     ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
81     ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
82     ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
83     ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
84     /* put in numbers */
85     base = (dd->base)/dd->w;
86     for (y=ymin; y<=ymax; y++) {
87       for (x=xmin; x<=xmax; x++) {
88         ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)base++);CHKERRQ(ierr);
89         ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
90       }
91     }
92     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
93 
94     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
95     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
96 
97     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
98     /* overlay ghost numbers, useful for error checking */
99     ierr = ISLocalToGlobalMappingGetBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr);
100     base = 0; xmin = dd->Xs; xmax = dd->Xe; ymin = dd->Ys; ymax = dd->Ye;
101     for (y=ymin; y<ymax; y++) {
102       for (x=xmin; x<xmax; x++) {
103         if ((base % dd->w) == 0) {
104           ierr = PetscSNPrintf(node,sizeof(node),"%d",(int)(idx[base/dd->w]));CHKERRQ(ierr);
105           ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
106         }
107         base++;
108       }
109     }
110     ierr = ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap,&idx);CHKERRQ(ierr);
111     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
112 
113     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
114     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
115   } else if (isbinary) {
116     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
117 #if defined(PETSC_HAVE_MATLAB_ENGINE)
118   } else if (ismatlab) {
119     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
120 #endif
121   }
122   PetscFunctionReturn(0);
123 }
124 
125 /*
126       M is number of grid points
127       m is number of processors
128 
129 */
130 #undef __FUNCT__
131 #define __FUNCT__ "DMDASplitComm2d"
132 PetscErrorCode  DMDASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm)
133 {
134   PetscErrorCode ierr;
135   PetscInt       m,n = 0,x = 0,y = 0;
136   PetscMPIInt    size,csize,rank;
137 
138   PetscFunctionBegin;
139   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
140   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
141 
142   csize = 4*size;
143   do {
144     if (csize % 4) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot split communicator of size %d tried %d %D %D",size,csize,x,y);
145     csize = csize/4;
146 
147     m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)csize)/((PetscReal)N)));
148     if (!m) m = 1;
149     while (m > 0) {
150       n = csize/m;
151       if (m*n == csize) break;
152       m--;
153     }
154     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
155 
156     x = M/m + ((M % m) > ((csize-1) % m));
157     y = (N + (csize-1)/m)/n;
158   } while ((x < 4 || y < 4) && csize > 1);
159   if (size != csize) {
160     MPI_Group   entire_group,sub_group;
161     PetscMPIInt i,*groupies;
162 
163     ierr = MPI_Comm_group(comm,&entire_group);CHKERRQ(ierr);
164     ierr = PetscMalloc1(csize,&groupies);CHKERRQ(ierr);
165     for (i=0; i<csize; i++) {
166       groupies[i] = (rank/csize)*csize + i;
167     }
168     ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);CHKERRQ(ierr);
169     ierr = PetscFree(groupies);CHKERRQ(ierr);
170     ierr = MPI_Comm_create(comm,sub_group,outcomm);CHKERRQ(ierr);
171     ierr = MPI_Group_free(&entire_group);CHKERRQ(ierr);
172     ierr = MPI_Group_free(&sub_group);CHKERRQ(ierr);
173     ierr = PetscInfo1(0,"DMDASplitComm2d:Creating redundant coarse problems of size %d\n",csize);CHKERRQ(ierr);
174   } else {
175     *outcomm = comm;
176   }
177   PetscFunctionReturn(0);
178 }
179 
180 #if defined(new)
181 #undef __FUNCT__
182 #define __FUNCT__ "DMDAGetDiagonal_MFFD"
183 /*
184   DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
185     function lives on a DMDA
186 
187         y ~= (F(u + ha) - F(u))/h,
188   where F = nonlinear function, as set by SNESSetFunction()
189         u = current iterate
190         h = difference interval
191 */
192 PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a)
193 {
194   PetscScalar    h,*aa,*ww,v;
195   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
196   PetscErrorCode ierr;
197   PetscInt       gI,nI;
198   MatStencil     stencil;
199   DMDALocalInfo  info;
200 
201   PetscFunctionBegin;
202   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
203   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
204 
205   ierr = VecGetArray(U,&ww);CHKERRQ(ierr);
206   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
207 
208   nI = 0;
209   h  = ww[gI];
210   if (h == 0.0) h = 1.0;
211   if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin;
212   else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
213   h *= epsilon;
214 
215   ww[gI] += h;
216   ierr    = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
217   aa[nI]  = (v - aa[nI])/h;
218   ww[gI] -= h;
219   nI++;
220 
221   ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr);
222   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
223   PetscFunctionReturn(0);
224 }
225 #endif
226 
227 #undef __FUNCT__
228 #define __FUNCT__ "DMSetUp_DA_2D"
229 PetscErrorCode  DMSetUp_DA_2D(DM da)
230 {
231   DM_DA            *dd = (DM_DA*)da->data;
232   const PetscInt   M            = dd->M;
233   const PetscInt   N            = dd->N;
234   PetscInt         m            = dd->m;
235   PetscInt         n            = dd->n;
236   const PetscInt   dof          = dd->w;
237   const PetscInt   s            = dd->s;
238   DMBoundaryType   bx           = dd->bx;
239   DMBoundaryType   by           = dd->by;
240   DMDAStencilType  stencil_type = dd->stencil_type;
241   PetscInt         *lx          = dd->lx;
242   PetscInt         *ly          = dd->ly;
243   MPI_Comm         comm;
244   PetscMPIInt      rank,size;
245   PetscInt         xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,IXs,IXe,IYs,IYe;
246   PetscInt         up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn;
247   PetscInt         xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
248   PetscInt         s_x,s_y; /* s proportionalized to w */
249   PetscInt         sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
250   Vec              local,global;
251   VecScatter       gtol;
252   IS               to,from;
253   PetscErrorCode   ierr;
254 
255   PetscFunctionBegin;
256   if (stencil_type == DMDA_STENCIL_BOX && (bx == DM_BOUNDARY_MIRROR || by == DM_BOUNDARY_MIRROR)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Mirror boundary and box stencil");
257   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
258 #if !defined(PETSC_USE_64BIT_INDICES)
259   if (((Petsc64bitInt) M)*((Petsc64bitInt) N)*((Petsc64bitInt) dof) > (Petsc64bitInt) PETSC_MPI_INT_MAX) SETERRQ3(comm,PETSC_ERR_INT_OVERFLOW,"Mesh of %D by %D by %D (dof) is too large for 32 bit indices",M,N,dof);
260 #endif
261 
262   if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
263   if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
264 
265   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
266   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
267 
268   dd->p = 1;
269   if (m != PETSC_DECIDE) {
270     if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
271     else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
272   }
273   if (n != PETSC_DECIDE) {
274     if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
275     else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
276   }
277 
278   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
279     if (n != PETSC_DECIDE) {
280       m = size/n;
281     } else if (m != PETSC_DECIDE) {
282       n = size/m;
283     } else {
284       /* try for squarish distribution */
285       m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N)));
286       if (!m) m = 1;
287       while (m > 0) {
288         n = size/m;
289         if (m*n == size) break;
290         m--;
291       }
292       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
293     }
294     if (m*n != size) SETERRQ(comm,PETSC_ERR_PLIB,"Unable to create partition, check the size of the communicator and input m and n ");
295   } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
296 
297   if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
298   if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
299 
300   /*
301      Determine locally owned region
302      xs is the first local node number, x is the number of local nodes
303   */
304   if (!lx) {
305     ierr = PetscMalloc1(m, &dd->lx);CHKERRQ(ierr);
306     lx   = dd->lx;
307     for (i=0; i<m; i++) {
308       lx[i] = M/m + ((M % m) > i);
309     }
310   }
311   x  = lx[rank % m];
312   xs = 0;
313   for (i=0; i<(rank % m); i++) {
314     xs += lx[i];
315   }
316 #if defined(PETSC_USE_DEBUG)
317   left = xs;
318   for (i=(rank % m); i<m; i++) {
319     left += lx[i];
320   }
321   if (left != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M: %D %D",left,M);
322 #endif
323 
324   /*
325      Determine locally owned region
326      ys is the first local node number, y is the number of local nodes
327   */
328   if (!ly) {
329     ierr = PetscMalloc1(n, &dd->ly);CHKERRQ(ierr);
330     ly   = dd->ly;
331     for (i=0; i<n; i++) {
332       ly[i] = N/n + ((N % n) > i);
333     }
334   }
335   y  = ly[rank/m];
336   ys = 0;
337   for (i=0; i<(rank/m); i++) {
338     ys += ly[i];
339   }
340 #if defined(PETSC_USE_DEBUG)
341   left = ys;
342   for (i=(rank/m); i<n; i++) {
343     left += ly[i];
344   }
345   if (left != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of ly across processors not equal to N: %D %D",left,N);
346 #endif
347 
348   /*
349    check if the scatter requires more than one process neighbor or wraps around
350    the domain more than once
351   */
352   if ((x < s) && ((m > 1) || (bx == DM_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s);
353   if ((y < s) && ((n > 1) || (by == DM_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s);
354   xe = xs + x;
355   ye = ys + y;
356 
357   /* determine ghost region (Xs) and region scattered into (IXs)  */
358   if (xs-s > 0) {
359     Xs = xs - s; IXs = xs - s;
360   } else {
361     if (bx) {
362       Xs = xs - s;
363     } else {
364       Xs = 0;
365     }
366     IXs = 0;
367   }
368   if (xe+s <= M) {
369     Xe = xe + s; IXe = xe + s;
370   } else {
371     if (bx) {
372       Xs = xs - s; Xe = xe + s;
373     } else {
374       Xe = M;
375     }
376     IXe = M;
377   }
378 
379   if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) {
380     IXs = xs - s;
381     IXe = xe + s;
382     Xs  = xs - s;
383     Xe  = xe + s;
384   }
385 
386   if (ys-s > 0) {
387     Ys = ys - s; IYs = ys - s;
388   } else {
389     if (by) {
390       Ys = ys - s;
391     } else {
392       Ys = 0;
393     }
394     IYs = 0;
395   }
396   if (ye+s <= N) {
397     Ye = ye + s; IYe = ye + s;
398   } else {
399     if (by) {
400       Ye = ye + s;
401     } else {
402       Ye = N;
403     }
404     IYe = N;
405   }
406 
407   if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) {
408     IYs = ys - s;
409     IYe = ye + s;
410     Ys  = ys - s;
411     Ye  = ye + s;
412   }
413 
414   /* stencil length in each direction */
415   s_x = s;
416   s_y = s;
417 
418   /* determine starting point of each processor */
419   nn       = x*y;
420   ierr     = PetscMalloc2(size+1,&bases,size,&ldims);CHKERRQ(ierr);
421   ierr     = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
422   bases[0] = 0;
423   for (i=1; i<=size; i++) {
424     bases[i] = ldims[i-1];
425   }
426   for (i=1; i<=size; i++) {
427     bases[i] += bases[i-1];
428   }
429   base = bases[rank]*dof;
430 
431   /* allocate the base parallel and sequential vectors */
432   dd->Nlocal = x*y*dof;
433   ierr       = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,NULL,&global);CHKERRQ(ierr);
434   dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof;
435   ierr       = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,NULL,&local);CHKERRQ(ierr);
436 
437   /* generate appropriate vector scatters */
438   /* local to global inserts non-ghost point region into global */
439   ierr  = PetscMalloc1((IXe-IXs)*(IYe-IYs),&idx);CHKERRQ(ierr);
440   left  = xs - Xs; right = left + x;
441   down  = ys - Ys; up = down + y;
442   count = 0;
443   for (i=down; i<up; i++) {
444     for (j=left; j<right; j++) {
445       idx[count++] = i*(Xe-Xs) + j;
446     }
447   }
448 
449   /* global to local must include ghost points within the domain,
450      but not ghost points outside the domain that aren't periodic */
451   if (stencil_type == DMDA_STENCIL_BOX) {
452     left  = IXs - Xs; right = left + (IXe-IXs);
453     down  = IYs - Ys; up = down + (IYe-IYs);
454     count = 0;
455     for (i=down; i<up; i++) {
456       for (j=left; j<right; j++) {
457         idx[count++] = j + i*(Xe-Xs);
458       }
459     }
460     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
461 
462   } else {
463     /* must drop into cross shape region */
464     /*       ---------|
465             |  top    |
466          |---         ---| up
467          |   middle      |
468          |               |
469          ----         ---- down
470             | bottom  |
471             -----------
472          Xs xs        xe Xe */
473     left  = xs - Xs; right = left + x;
474     down  = ys - Ys; up = down + y;
475     count = 0;
476     /* bottom */
477     for (i=(IYs-Ys); i<down; i++) {
478       for (j=left; j<right; j++) {
479         idx[count++] = j + i*(Xe-Xs);
480       }
481     }
482     /* middle */
483     for (i=down; i<up; i++) {
484       for (j=(IXs-Xs); j<(IXe-Xs); j++) {
485         idx[count++] = j + i*(Xe-Xs);
486       }
487     }
488     /* top */
489     for (i=up; i<up+IYe-ye; i++) {
490       for (j=left; j<right; j++) {
491         idx[count++] = j + i*(Xe-Xs);
492       }
493     }
494     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
495   }
496 
497 
498   /* determine who lies on each side of us stored in    n6 n7 n8
499                                                         n3    n5
500                                                         n0 n1 n2
501   */
502 
503   /* Assume the Non-Periodic Case */
504   n1 = rank - m;
505   if (rank % m) {
506     n0 = n1 - 1;
507   } else {
508     n0 = -1;
509   }
510   if ((rank+1) % m) {
511     n2 = n1 + 1;
512     n5 = rank + 1;
513     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
514   } else {
515     n2 = -1; n5 = -1; n8 = -1;
516   }
517   if (rank % m) {
518     n3 = rank - 1;
519     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
520   } else {
521     n3 = -1; n6 = -1;
522   }
523   n7 = rank + m; if (n7 >= m*n) n7 = -1;
524 
525   if (bx == DM_BOUNDARY_PERIODIC && by == DM_BOUNDARY_PERIODIC) {
526     /* Modify for Periodic Cases */
527     /* Handle all four corners */
528     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
529     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
530     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
531     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
532 
533     /* Handle Top and Bottom Sides */
534     if (n1 < 0) n1 = rank + m * (n-1);
535     if (n7 < 0) n7 = rank - m * (n-1);
536     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
537     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
538     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
539     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
540 
541     /* Handle Left and Right Sides */
542     if (n3 < 0) n3 = rank + (m-1);
543     if (n5 < 0) n5 = rank - (m-1);
544     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
545     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
546     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
547     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
548   } else if (by == DM_BOUNDARY_PERIODIC) {  /* Handle Top and Bottom Sides */
549     if (n1 < 0) n1 = rank + m * (n-1);
550     if (n7 < 0) n7 = rank - m * (n-1);
551     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
552     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
553     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
554     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
555   } else if (bx == DM_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */
556     if (n3 < 0) n3 = rank + (m-1);
557     if (n5 < 0) n5 = rank - (m-1);
558     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
559     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
560     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
561     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
562   }
563 
564   ierr = PetscMalloc1(9,&dd->neighbors);CHKERRQ(ierr);
565 
566   dd->neighbors[0] = n0;
567   dd->neighbors[1] = n1;
568   dd->neighbors[2] = n2;
569   dd->neighbors[3] = n3;
570   dd->neighbors[4] = rank;
571   dd->neighbors[5] = n5;
572   dd->neighbors[6] = n6;
573   dd->neighbors[7] = n7;
574   dd->neighbors[8] = n8;
575 
576   if (stencil_type == DMDA_STENCIL_STAR) {
577     /* save corner processor numbers */
578     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
579     n0  = n2 = n6 = n8 = -1;
580   }
581 
582   ierr = PetscMalloc1((Xe-Xs)*(Ye-Ys),&idx);CHKERRQ(ierr);
583 
584   nn = 0;
585   xbase = bases[rank];
586   for (i=1; i<=s_y; i++) {
587     if (n0 >= 0) { /* left below */
588       x_t = lx[n0 % m];
589       y_t = ly[(n0/m)];
590       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
591       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
592     }
593 
594     if (n1 >= 0) { /* directly below */
595       x_t = x;
596       y_t = ly[(n1/m)];
597       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
598       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
599     } else if (by == DM_BOUNDARY_MIRROR) {
600       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
601     }
602 
603     if (n2 >= 0) { /* right below */
604       x_t = lx[n2 % m];
605       y_t = ly[(n2/m)];
606       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
607       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
608     }
609   }
610 
611   for (i=0; i<y; i++) {
612     if (n3 >= 0) { /* directly left */
613       x_t = lx[n3 % m];
614       /* y_t = y; */
615       s_t = bases[n3] + (i+1)*x_t - s_x;
616       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
617     } else if (bx == DM_BOUNDARY_MIRROR) {
618       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
619     }
620 
621     for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
622 
623     if (n5 >= 0) { /* directly right */
624       x_t = lx[n5 % m];
625       /* y_t = y; */
626       s_t = bases[n5] + (i)*x_t;
627       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
628     } else if (bx == DM_BOUNDARY_MIRROR) {
629       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
630     }
631   }
632 
633   for (i=1; i<=s_y; i++) {
634     if (n6 >= 0) { /* left above */
635       x_t = lx[n6 % m];
636       /* y_t = ly[(n6/m)]; */
637       s_t = bases[n6] + (i)*x_t - s_x;
638       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
639     }
640 
641     if (n7 >= 0) { /* directly above */
642       x_t = x;
643       /* y_t = ly[(n7/m)]; */
644       s_t = bases[n7] + (i-1)*x_t;
645       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
646     } else if (by == DM_BOUNDARY_MIRROR) {
647       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
648     }
649 
650     if (n8 >= 0) { /* right above */
651       x_t = lx[n8 % m];
652       /* y_t = ly[(n8/m)]; */
653       s_t = bases[n8] + (i-1)*x_t;
654       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
655     }
656   }
657 
658   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_USE_POINTER,&from);CHKERRQ(ierr);
659   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
660   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);CHKERRQ(ierr);
661   ierr = ISDestroy(&to);CHKERRQ(ierr);
662   ierr = ISDestroy(&from);CHKERRQ(ierr);
663 
664   if (stencil_type == DMDA_STENCIL_STAR) {
665     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
666   }
667 
668   if (((stencil_type == DMDA_STENCIL_STAR)  ||
669        (bx && bx != DM_BOUNDARY_PERIODIC) ||
670        (by && by != DM_BOUNDARY_PERIODIC))) {
671     /*
672         Recompute the local to global mappings, this time keeping the
673       information about the cross corner processor numbers and any ghosted
674       but not periodic indices.
675     */
676     nn    = 0;
677     xbase = bases[rank];
678     for (i=1; i<=s_y; i++) {
679       if (n0 >= 0) { /* left below */
680         x_t = lx[n0 % m];
681         y_t = ly[(n0/m)];
682         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
683         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
684       } else if (xs-Xs > 0 && ys-Ys > 0) {
685         for (j=0; j<s_x; j++) idx[nn++] = -1;
686       }
687       if (n1 >= 0) { /* directly below */
688         x_t = x;
689         y_t = ly[(n1/m)];
690         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
691         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
692       } else if (ys-Ys > 0) {
693         if (by == DM_BOUNDARY_MIRROR) {
694           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
695         } else {
696           for (j=0; j<x; j++) idx[nn++] = -1;
697         }
698       }
699       if (n2 >= 0) { /* right below */
700         x_t = lx[n2 % m];
701         y_t = ly[(n2/m)];
702         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
703         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
704       } else if (Xe-xe> 0 && ys-Ys > 0) {
705         for (j=0; j<s_x; j++) idx[nn++] = -1;
706       }
707     }
708 
709     for (i=0; i<y; i++) {
710       if (n3 >= 0) { /* directly left */
711         x_t = lx[n3 % m];
712         /* y_t = y; */
713         s_t = bases[n3] + (i+1)*x_t - s_x;
714         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
715       } else if (xs-Xs > 0) {
716         if (bx == DM_BOUNDARY_MIRROR) {
717           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
718         } else {
719           for (j=0; j<s_x; j++) idx[nn++] = -1;
720         }
721       }
722 
723       for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */
724 
725       if (n5 >= 0) { /* directly right */
726         x_t = lx[n5 % m];
727         /* y_t = y; */
728         s_t = bases[n5] + (i)*x_t;
729         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
730       } else if (Xe-xe > 0) {
731         if (bx == DM_BOUNDARY_MIRROR) {
732           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
733         } else {
734           for (j=0; j<s_x; j++) idx[nn++] = -1;
735         }
736       }
737     }
738 
739     for (i=1; i<=s_y; i++) {
740       if (n6 >= 0) { /* left above */
741         x_t = lx[n6 % m];
742         /* y_t = ly[(n6/m)]; */
743         s_t = bases[n6] + (i)*x_t - s_x;
744         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
745       } else if (xs-Xs > 0 && Ye-ye > 0) {
746         for (j=0; j<s_x; j++) idx[nn++] = -1;
747       }
748       if (n7 >= 0) { /* directly above */
749         x_t = x;
750         /* y_t = ly[(n7/m)]; */
751         s_t = bases[n7] + (i-1)*x_t;
752         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
753       } else if (Ye-ye > 0) {
754         if (by == DM_BOUNDARY_MIRROR) {
755           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
756         } else {
757           for (j=0; j<x; j++) idx[nn++] = -1;
758         }
759       }
760       if (n8 >= 0) { /* right above */
761         x_t = lx[n8 % m];
762         /* y_t = ly[(n8/m)]; */
763         s_t = bases[n8] + (i-1)*x_t;
764         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
765       } else if (Xe-xe > 0 && Ye-ye > 0) {
766         for (j=0; j<s_x; j++) idx[nn++] = -1;
767       }
768     }
769   }
770   /*
771      Set the local to global ordering in the global vector, this allows use
772      of VecSetValuesLocal().
773   */
774   ierr = ISLocalToGlobalMappingCreate(comm,dof,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr);
775   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);CHKERRQ(ierr);
776 
777   ierr  = PetscFree2(bases,ldims);CHKERRQ(ierr);
778   dd->m = m;  dd->n  = n;
779   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
780   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
781   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
782 
783   ierr = VecDestroy(&local);CHKERRQ(ierr);
784   ierr = VecDestroy(&global);CHKERRQ(ierr);
785 
786   dd->gtol      = gtol;
787   dd->base      = base;
788   da->ops->view = DMView_DA_2d;
789   dd->ltol      = NULL;
790   dd->ao        = NULL;
791   PetscFunctionReturn(0);
792 }
793 
794 #undef __FUNCT__
795 #define __FUNCT__ "DMDACreate2d"
796 /*@C
797    DMDACreate2d -  Creates an object that will manage the communication of  two-dimensional
798    regular array data that is distributed across some processors.
799 
800    Collective on MPI_Comm
801 
802    Input Parameters:
803 +  comm - MPI communicator
804 .  bx,by - type of ghost nodes the array have.
805          Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC.
806 .  stencil_type - stencil type.  Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
807 .  M,N - global dimension in each direction of the array (use -M and or -N to indicate that it may be set to a different value
808             from the command line with -da_grid_x <M> -da_grid_y <N>)
809 .  m,n - corresponding number of processors in each dimension
810          (or PETSC_DECIDE to have calculated)
811 .  dof - number of degrees of freedom per node
812 .  s - stencil width
813 -  lx, ly - arrays containing the number of nodes in each cell along
814            the x and y coordinates, or NULL. If non-null, these
815            must be of length as m and n, and the corresponding
816            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
817            must be M, and the sum of the ly[] entries must be N.
818 
819    Output Parameter:
820 .  da - the resulting distributed array object
821 
822    Options Database Key:
823 +  -dm_view - Calls DMView() at the conclusion of DMDACreate2d()
824 .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
825 .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
826 .  -da_processors_x <nx> - number of processors in x direction
827 .  -da_processors_y <ny> - number of processors in y direction
828 .  -da_refine_x <rx> - refinement ratio in x direction
829 .  -da_refine_y <ry> - refinement ratio in y direction
830 -  -da_refine <n> - refine the DMDA n times before creating, if M or N < 0
831 
832 
833    Level: beginner
834 
835    Notes:
836    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
837    standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
838    the standard 9-pt stencil.
839 
840    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
841    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
842    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
843 
844 .keywords: distributed array, create, two-dimensional
845 
846 .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
847           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(),
848           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
849 
850 @*/
851 
852 PetscErrorCode  DMDACreate2d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMDAStencilType stencil_type,
853                           PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
854 {
855   PetscErrorCode ierr;
856 
857   PetscFunctionBegin;
858   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
859   ierr = DMSetDimension(*da, 2);CHKERRQ(ierr);
860   ierr = DMDASetSizes(*da, M, N, 1);CHKERRQ(ierr);
861   ierr = DMDASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr);
862   ierr = DMDASetBoundaryType(*da, bx, by, DM_BOUNDARY_NONE);CHKERRQ(ierr);
863   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
864   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
865   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
866   ierr = DMDASetOwnershipRanges(*da, lx, ly, NULL);CHKERRQ(ierr);
867   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
868   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
869   ierr = DMSetUp(*da);CHKERRQ(ierr);
870   PetscFunctionReturn(0);
871 }
872