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