xref: /petsc/src/dm/impls/da/da2.c (revision fa59e80594ae4f5b96fb1a0d88f82abc4bbaec69)
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 !defined(PETSC_USE_COMPLEX)
202     if (h < umin && h >= 0.0)      h = umin;
203     else if (h < 0.0 && h > -umin) h = -umin;
204 #else
205     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
206     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
207 #endif
208     h     *= epsilon;
209 
210     ww[gI] += h;
211     ierr          = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
212     aa[nI]  = (v - aa[nI])/h;
213     ww[gI] -= h;
214     nI++;
215   }
216   ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr);
217   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
218   PetscFunctionReturn(0);
219 }
220 #endif
221 
222 #undef __FUNCT__
223 #define __FUNCT__ "DMSetUp_DA_2D"
224 PetscErrorCode  DMSetUp_DA_2D(DM da)
225 {
226   DM_DA            *dd = (DM_DA*)da->data;
227   const PetscInt   M            = dd->M;
228   const PetscInt   N            = dd->N;
229   PetscInt         m            = dd->m;
230   PetscInt         n            = dd->n;
231   PetscInt         o            = dd->overlap;
232   const PetscInt   dof          = dd->w;
233   const PetscInt   s            = dd->s;
234   DMDABoundaryType bx         = dd->bx;
235   DMDABoundaryType by         = dd->by;
236   DMDAStencilType  stencil_type = dd->stencil_type;
237   PetscInt         *lx           = dd->lx;
238   PetscInt         *ly           = dd->ly;
239   MPI_Comm         comm;
240   PetscMPIInt      rank,size;
241   PetscInt         xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,start,end,IXs,IXe,IYs,IYe;
242   PetscInt         up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn,*idx_cpy;
243   const PetscInt   *idx_full;
244   PetscInt         xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
245   PetscInt         s_x,s_y; /* s proportionalized to w */
246   PetscInt         sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
247   Vec              local,global;
248   VecScatter       ltog,gtol;
249   IS               to,from,ltogis;
250   PetscErrorCode   ierr;
251 
252   PetscFunctionBegin;
253   if (dof < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
254   if (s < 0) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
255   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
256 #if !defined(PETSC_USE_64BIT_INDICES)
257   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);
258 #endif
259 
260   if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
261   if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
262 
263   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
264   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
265 
266   dd->dim         = 2;
267   ierr = PetscMalloc(dof*sizeof(char*),&dd->fieldname);CHKERRQ(ierr);
268   ierr = PetscMemzero(dd->fieldname,dof*sizeof(char*));CHKERRQ(ierr);
269 
270   if (m != PETSC_DECIDE) {
271     if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
272     else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
273   }
274   if (n != PETSC_DECIDE) {
275     if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
276     else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
277   }
278 
279   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
280     if (n != PETSC_DECIDE) {
281       m = size/n;
282     } else if (m != PETSC_DECIDE) {
283       n = size/m;
284     } else {
285       /* try for squarish distribution */
286       m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)N)));
287       if (!m) m = 1;
288       while (m > 0) {
289 	n = size/m;
290 	if (m*n == size) break;
291 	m--;
292       }
293       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
294     }
295     if (m*n != size) SETERRQ(comm,PETSC_ERR_PLIB,"Unable to create partition, check the size of the communicator and input m and n ");
296   } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
297 
298   if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
299   if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
300 
301   /*
302      Determine locally owned region
303      xs is the first local node number, x is the number of local nodes
304   */
305   if (!lx) {
306     ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
307     lx = dd->lx;
308     for (i=0; i<m; i++) {
309       lx[i] = M/m + ((M % m) > i);
310     }
311   }
312   x  = lx[rank % m];
313   xs = 0;
314   for (i=0; i<(rank % m); i++) {
315     xs += lx[i];
316   }
317 #if defined(PETSC_USE_DEBUG)
318   left = xs;
319   for (i=(rank % m); i<m; i++) {
320     left += lx[i];
321   }
322   if (left != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M: %D %D",left,M);
323 #endif
324 
325   /*
326      Determine locally owned region
327      ys is the first local node number, y is the number of local nodes
328   */
329   if (!ly) {
330     ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
331     ly = dd->ly;
332     for (i=0; i<n; i++) {
333       ly[i] = N/n + ((N % n) > i);
334     }
335   }
336   y  = ly[rank/m];
337   ys = 0;
338   for (i=0; i<(rank/m); i++) {
339     ys += ly[i];
340   }
341 #if defined(PETSC_USE_DEBUG)
342   left = ys;
343   for (i=(rank/m); i<n; i++) {
344     left += ly[i];
345   }
346   if (left != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of ly across processors not equal to N: %D %D",left,N);
347 #endif
348 
349   /*
350    check if the scatter requires more than one process neighbor or wraps around
351    the domain more than once
352   */
353   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);
354   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);
355   xe = xs + x;
356   ye = ys + y;
357 
358   /* determine ghost region (Xs) and region scattered into (IXs)  */
359   if (xs-s-o > 0) {
360     Xs = xs - s - o; IXs = xs - s - o;
361   } else {
362     if (bx) {
363       Xs = xs - s;
364     } else {
365       Xs = 0;
366     }
367     IXs = 0;
368   }
369   if (xe+s+o <= M) {
370     Xe = xe + s + o; IXe = xe + s + o;
371   } else {
372     if (bx) {
373       Xs = xs - s - o; Xe = xe + s;
374     } else {
375       Xe = M;
376     }
377     IXe = M;
378   }
379 
380   if (bx == DMDA_BOUNDARY_PERIODIC) {
381     IXs = xs - s - o;
382     IXe = xe + s + o;
383     Xs = xs - s - o;
384     Xe = xe + s + o;
385   }
386 
387   if (ys-s-o > 0) {
388     Ys = ys - s - o; IYs = ys - s - o;
389   } else {
390     if (by) {
391       Ys = ys - s;
392     } else {
393       Ys = 0;
394     }
395     IYs = 0;
396   }
397   if (ye+s+o <= N) {
398     Ye = ye + s + o; IYe = ye + s + o;
399   } else {
400     if (by) {
401       Ye = ye + s;
402     } else {
403       Ye = N;
404     }
405     IYe = N;
406   }
407 
408   if (by == DMDA_BOUNDARY_PERIODIC) {
409     IYs = ys - s - o;
410     IYe = ye + s + o;
411     Ys = ys - s - o;
412     Ye = ye + s + o;
413   }
414 
415   /* stencil length in each direction */
416   s_x = s + o;
417   s_y = s + o;
418 
419   /* determine starting point of each processor */
420   nn    = x*y;
421   ierr  = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr);
422   ierr  = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
423   bases[0] = 0;
424   for (i=1; i<=size; i++) {
425     bases[i] = ldims[i-1];
426   }
427   for (i=1; i<=size; i++) {
428     bases[i] += bases[i-1];
429   }
430   base = bases[rank]*dof;
431 
432   /* allocate the base parallel and sequential vectors */
433   dd->Nlocal = x*y*dof;
434   ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr);
435   dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof;
436   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,0,&local);CHKERRQ(ierr);
437 
438   /* generate appropriate vector scatters */
439   /* local to global inserts non-ghost point region into global */
440   ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr);
441   ierr = ISCreateStride(comm,x*y*dof,start,1,&to);CHKERRQ(ierr);
442 
443   count = x*y;
444   ierr = PetscMalloc(x*y*sizeof(PetscInt),&idx);CHKERRQ(ierr);
445   left = xs - Xs; right = left + x;
446   down = ys - Ys; up = down + y;
447   count = 0;
448   for (i=down; i<up; i++) {
449     for (j=left; j<right; j++) {
450       idx[count++] = i*(Xe-Xs) + j;
451     }
452   }
453 
454   ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
455   ierr = VecScatterCreate(local,from,global,to,&ltog);CHKERRQ(ierr);
456   ierr = PetscLogObjectParent(dd,ltog);CHKERRQ(ierr);
457   ierr = ISDestroy(&from);CHKERRQ(ierr);
458   ierr = ISDestroy(&to);CHKERRQ(ierr);
459 
460   /* global to local must include ghost points within the domain,
461      but not ghost points outside the domain that aren't periodic */
462   if (stencil_type == DMDA_STENCIL_BOX || o > 0) {
463     count = (IXe-IXs)*(IYe-IYs);
464     ierr  = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
465 
466     left = IXs - Xs; right = left + (IXe-IXs);
467     down = IYs - Ys; up = down + (IYe-IYs);
468     count = 0;
469     for (i=down; i<up; i++) {
470       for (j=left; j<right; j++) {
471         idx[count++] = j + i*(Xe-Xs);
472       }
473     }
474     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
475 
476   } else {
477     /* must drop into cross shape region */
478     /*       ---------|
479             |  top    |
480          |---         ---| up
481          |   middle      |
482          |               |
483          ----         ---- down
484             | bottom  |
485             -----------
486          Xs xs        xe Xe */
487     count = (ys-IYs)*x + y*(IXe-IXs) + (IYe-ye)*x;
488     ierr  = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
489 
490     left = xs - Xs; right = left + x;
491     down = ys - Ys; up = down + y;
492     count = 0;
493     /* bottom */
494     for (i=(IYs-Ys); i<down; i++) {
495       for (j=left; j<right; j++) {
496         idx[count++] = j + i*(Xe-Xs);
497       }
498     }
499     /* middle */
500     for (i=down; i<up; i++) {
501       for (j=(IXs-Xs); j<(IXe-Xs); j++) {
502         idx[count++] = j + i*(Xe-Xs);
503       }
504     }
505     /* top */
506     for (i=up; i<up+IYe-ye; i++) {
507       for (j=left; j<right; j++) {
508         idx[count++] = j + i*(Xe-Xs);
509       }
510     }
511     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
512   }
513 
514 
515   /* determine who lies on each side of us stored in    n6 n7 n8
516                                                         n3    n5
517                                                         n0 n1 n2
518   */
519 
520   /* Assume the Non-Periodic Case */
521   n1 = rank - m;
522   if (rank % m) {
523     n0 = n1 - 1;
524   } else {
525     n0 = -1;
526   }
527   if ((rank+1) % m) {
528     n2 = n1 + 1;
529     n5 = rank + 1;
530     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
531   } else {
532     n2 = -1; n5 = -1; n8 = -1;
533   }
534   if (rank % m) {
535     n3 = rank - 1;
536     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
537   } else {
538     n3 = -1; n6 = -1;
539   }
540   n7 = rank + m; if (n7 >= m*n) n7 = -1;
541 
542   if (bx == DMDA_BOUNDARY_PERIODIC && by == DMDA_BOUNDARY_PERIODIC) {
543   /* Modify for Periodic Cases */
544     /* Handle all four corners */
545     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
546     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
547     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
548     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
549 
550     /* 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 
558     /* Handle Left and Right Sides */
559     if (n3 < 0) n3 = rank + (m-1);
560     if (n5 < 0) n5 = rank - (m-1);
561     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
562     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
563     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
564     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
565   } else if (by == DMDA_BOUNDARY_PERIODIC) {  /* Handle Top and Bottom Sides */
566     if (n1 < 0) n1 = rank + m * (n-1);
567     if (n7 < 0) n7 = rank - m * (n-1);
568     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
569     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
570     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
571     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
572   } else if (bx == DMDA_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */
573     if (n3 < 0) n3 = rank + (m-1);
574     if (n5 < 0) n5 = rank - (m-1);
575     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
576     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
577     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
578     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
579   }
580 
581   ierr = PetscMalloc(9*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr);
582   dd->neighbors[0] = n0;
583   dd->neighbors[1] = n1;
584   dd->neighbors[2] = n2;
585   dd->neighbors[3] = n3;
586   dd->neighbors[4] = rank;
587   dd->neighbors[5] = n5;
588   dd->neighbors[6] = n6;
589   dd->neighbors[7] = n7;
590   dd->neighbors[8] = n8;
591 
592   if (stencil_type == DMDA_STENCIL_STAR && o == 0) {
593     /* save corner processor numbers */
594     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
595     n0 = n2 = n6 = n8 = -1;
596   }
597 
598   ierr = PetscMalloc((Xe-Xs)*(Ye-Ys)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
599   ierr = PetscLogObjectMemory(da,(Xe-Xs)*(Ye-Ys)*sizeof(PetscInt));CHKERRQ(ierr);
600 
601   nn = 0;
602   xbase = bases[rank];
603   for (i=1; i<=s_y; i++) {
604     if (n0 >= 0) { /* left below */
605       x_t = lx[n0 % m];
606       y_t = ly[(n0/m)];
607       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
608       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
609     }
610     if (n1 >= 0) { /* directly below */
611       x_t = x;
612       y_t = ly[(n1/m)];
613       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
614       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
615     }
616     if (n2 >= 0) { /* right below */
617       x_t = lx[n2 % m];
618       y_t = ly[(n2/m)];
619       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
620       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
621     }
622   }
623 
624   for (i=0; i<y; i++) {
625     if (n3 >= 0) { /* directly left */
626       x_t = lx[n3 % m];
627       /* y_t = y; */
628       s_t = bases[n3] + (i+1)*x_t - s_x;
629       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
630     }
631 
632     for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
633 
634     if (n5 >= 0) { /* directly right */
635       x_t = lx[n5 % m];
636       /* y_t = y; */
637       s_t = bases[n5] + (i)*x_t;
638       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
639     }
640   }
641 
642   for (i=1; i<=s_y; i++) {
643     if (n6 >= 0) { /* left above */
644       x_t = lx[n6 % m];
645       /* y_t = ly[(n6/m)]; */
646       s_t = bases[n6] + (i)*x_t - s_x;
647       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
648     }
649     if (n7 >= 0) { /* directly above */
650       x_t = x;
651       /* y_t = ly[(n7/m)]; */
652       s_t = bases[n7] + (i-1)*x_t;
653       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
654     }
655     if (n8 >= 0) { /* right above */
656       x_t = lx[n8 % m];
657       /* y_t = ly[(n8/m)]; */
658       s_t = bases[n8] + (i-1)*x_t;
659       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
660     }
661   }
662 
663   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
664   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
665   ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr);
666   ierr = ISDestroy(&to);CHKERRQ(ierr);
667   ierr = ISDestroy(&from);CHKERRQ(ierr);
668 
669   if (stencil_type == DMDA_STENCIL_STAR && o == 0) {
670     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
671   }
672 
673   if (((stencil_type == DMDA_STENCIL_STAR) ||
674        (bx && bx != DMDA_BOUNDARY_PERIODIC) ||
675        (by && by != DMDA_BOUNDARY_PERIODIC)) && o == 0) {
676     /*
677         Recompute the local to global mappings, this time keeping the
678       information about the cross corner processor numbers and any ghosted
679       but not periodic indices.
680     */
681     nn = 0;
682     xbase = bases[rank];
683     for (i=1; i<=s_y; i++) {
684       if (n0 >= 0) { /* left below */
685         x_t = lx[n0 % m];
686         y_t = ly[(n0/m)];
687         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
688         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
689       } else if (xs-Xs > 0 && ys-Ys > 0) {
690         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
691       }
692       if (n1 >= 0) { /* directly below */
693         x_t = x;
694         y_t = ly[(n1/m)];
695         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
696         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
697       } else if (ys-Ys > 0) {
698         for (j=0; j<x; j++) { idx[nn++] = -1;}
699       }
700       if (n2 >= 0) { /* right below */
701         x_t = lx[n2 % m];
702         y_t = ly[(n2/m)];
703         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
704         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
705       } else if (Xe-xe> 0 && ys-Ys > 0) {
706         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
707       }
708     }
709 
710     for (i=0; i<y; i++) {
711       if (n3 >= 0) { /* directly left */
712         x_t = lx[n3 % m];
713         /* y_t = y; */
714         s_t = bases[n3] + (i+1)*x_t - s_x;
715         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
716       } else if (xs-Xs > 0) {
717         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
718       }
719 
720       for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
721 
722       if (n5 >= 0) { /* directly right */
723         x_t = lx[n5 % m];
724         /* y_t = y; */
725         s_t = bases[n5] + (i)*x_t;
726         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
727       } else if (Xe-xe > 0) {
728         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
729       }
730     }
731 
732     for (i=1; i<=s_y; i++) {
733       if (n6 >= 0) { /* left above */
734         x_t = lx[n6 % m];
735         /* y_t = ly[(n6/m)]; */
736         s_t = bases[n6] + (i)*x_t - s_x;
737         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
738       } else if (xs-Xs > 0 && Ye-ye > 0) {
739         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
740       }
741       if (n7 >= 0) { /* directly above */
742         x_t = x;
743         /* y_t = ly[(n7/m)]; */
744         s_t = bases[n7] + (i-1)*x_t;
745         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
746       } else if (Ye-ye > 0) {
747         for (j=0; j<x; j++) { idx[nn++] = -1;}
748       }
749       if (n8 >= 0) { /* right above */
750         x_t = lx[n8 % m];
751         /* y_t = ly[(n8/m)]; */
752         s_t = bases[n8] + (i-1)*x_t;
753         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
754       } else if (Xe-xe > 0 && Ye-ye > 0) {
755         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
756       }
757     }
758   }
759   /*
760      Set the local to global ordering in the global vector, this allows use
761      of VecSetValuesLocal().
762   */
763   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_OWN_POINTER,&ltogis);CHKERRQ(ierr);
764   ierr = PetscMalloc(nn*dof*sizeof(PetscInt),&idx_cpy);CHKERRQ(ierr);
765   ierr = PetscLogObjectMemory(da,nn*dof*sizeof(PetscInt));CHKERRQ(ierr);
766   ierr = ISGetIndices(ltogis, &idx_full);
767   ierr = PetscMemcpy(idx_cpy,idx_full,nn*dof*sizeof(PetscInt));CHKERRQ(ierr);
768   ierr = ISRestoreIndices(ltogis, &idx_full);
769   ierr = ISLocalToGlobalMappingCreateIS(ltogis,&da->ltogmap);CHKERRQ(ierr);
770   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
771   ierr = ISDestroy(&ltogis);CHKERRQ(ierr);
772   ierr = ISLocalToGlobalMappingBlock(da->ltogmap,dd->w,&da->ltogmapb);CHKERRQ(ierr);
773   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
774 
775   ierr = PetscFree2(bases,ldims);CHKERRQ(ierr);
776   dd->m  = m;  dd->n  = n;
777   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
778   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
779   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
780 
781   ierr = VecDestroy(&local);CHKERRQ(ierr);
782   ierr = VecDestroy(&global);CHKERRQ(ierr);
783 
784   dd->gtol      = gtol;
785   dd->ltog      = ltog;
786   dd->idx       = idx_cpy;
787   dd->Nl        = nn*dof;
788   dd->base      = base;
789   da->ops->view = DMView_DA_2d;
790   dd->ltol = PETSC_NULL;
791   dd->ao   = PETSC_NULL;
792 
793   PetscFunctionReturn(0);
794 }
795 
796 #undef __FUNCT__
797 #define __FUNCT__ "DMDACreate2d"
798 /*@C
799    DMDACreate2d -  Creates an object that will manage the communication of  two-dimensional
800    regular array data that is distributed across some processors.
801 
802    Collective on MPI_Comm
803 
804    Input Parameters:
805 +  comm - MPI communicator
806 .  bx,by - type of ghost nodes the array have.
807          Use one of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC.
808 .  stencil_type - stencil type.  Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
809 .  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
810             from the command line with -da_grid_x <M> -da_grid_y <N>)
811 .  m,n - corresponding number of processors in each dimension
812          (or PETSC_DECIDE to have calculated)
813 .  dof - number of degrees of freedom per node
814 .  s - stencil width
815 -  lx, ly - arrays containing the number of nodes in each cell along
816            the x and y coordinates, or PETSC_NULL. If non-null, these
817            must be of length as m and n, and the corresponding
818            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
819            must be M, and the sum of the ly[] entries must be N.
820 
821    Output Parameter:
822 .  da - the resulting distributed array object
823 
824    Options Database Key:
825 +  -da_view - Calls DMView() at the conclusion of DMDACreate2d()
826 .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
827 .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
828 .  -da_processors_x <nx> - number of processors in x direction
829 .  -da_processors_y <ny> - number of processors in y direction
830 .  -da_refine_x <rx> - refinement ratio in x direction
831 .  -da_refine_y <ry> - refinement ratio in y direction
832 -  -da_refine <n> - refine the DMDA n times before creating, if M or N < 0
833 
834 
835    Level: beginner
836 
837    Notes:
838    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
839    standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
840    the standard 9-pt stencil.
841 
842    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
843    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
844    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
845 
846 .keywords: distributed array, create, two-dimensional
847 
848 .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
849           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDASetRefinementFactor(),
850           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
851 
852 @*/
853 
854 PetscErrorCode  DMDACreate2d(MPI_Comm comm,DMDABoundaryType bx,DMDABoundaryType by,DMDAStencilType stencil_type,
855                           PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
856 {
857   PetscErrorCode ierr;
858 
859   PetscFunctionBegin;
860   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
861   ierr = DMDASetDim(*da, 2);CHKERRQ(ierr);
862   ierr = DMDASetSizes(*da, M, N, 1);CHKERRQ(ierr);
863   ierr = DMDASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr);
864   ierr = DMDASetBoundaryType(*da, bx, by, DMDA_BOUNDARY_NONE);CHKERRQ(ierr);
865   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
866   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
867   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
868   ierr = DMDASetOwnershipRanges(*da, lx, ly, PETSC_NULL);CHKERRQ(ierr);
869   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
870   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
871   ierr = DMSetUp(*da);CHKERRQ(ierr);
872   ierr = DMView_DA_Private(*da);CHKERRQ(ierr);
873   PetscFunctionReturn(0);
874 }
875