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