Actual source code: da2.c

petsc-dev 2014-02-02
Report Typos and Errors
  2: #include <petsc-private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
  3: #include <petscdraw.h>

  7: PetscErrorCode DMView_DA_2d(DM da,PetscViewer viewer)
  8: {
 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

 18:   MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);

 20:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
 21:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
 22:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
 23: #if defined(PETSC_HAVE_MATLAB_ENGINE)
 24:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);
 25: #endif
 26:   if (iascii) {
 27:     PetscViewerFormat format;

 29:     PetscViewerGetFormat(viewer, &format);
 30:     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) {
 31:       DMDALocalInfo info;
 32:       DMDAGetLocalInfo(da,&info);
 33:       PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);
 34:       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);
 35:       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);
 36:       PetscViewerFlush(viewer);
 37:       PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);
 38:     } else {
 39:       DMView_DA_VTK(da,viewer);
 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;

 51:     PetscViewerDrawGetDraw(viewer,0,&draw);
 52:     PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
 53:     if (!da->coordinates) {
 54:       PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);
 55:     }
 56:     PetscDrawSynchronizedClear(draw);

 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:         PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);
 63:       }
 64:       xmin = 0.0; xmax = dd->M - 1;
 65:       for (ymin=0; ymin<dd->N; ymin++) {
 66:         PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);
 67:       }
 68:     }
 69:     PetscDrawSynchronizedFlush(draw);
 70:     PetscDrawPause(draw);

 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:     PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);
 76:     PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);
 77:     PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);
 78:     PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);

 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:         PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);
 86:       }
 87:     }

 89:     PetscDrawSynchronizedFlush(draw);
 90:     PetscDrawPause(draw);
 91:     /* overlay ghost numbers, useful for error checking */
 92:     /* put in numbers */

 94:     base = 0;
 95:     ISLocalToGlobalMappingGetIndices(da->ltogmapb,&idx);
 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]));
101:           PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);
102:         }
103:         base++;
104:       }
105:     }
106:     ISLocalToGlobalMappingRestoreIndices(da->ltogmapb,&idx);
107:     PetscDrawSynchronizedFlush(draw);
108:     PetscDrawPause(draw);
109:   } else if (isbinary) {
110:     DMView_DA_Binary(da,viewer);
111: #if defined(PETSC_HAVE_MATLAB_ENGINE)
112:   } else if (ismatlab) {
113:     DMView_DA_Matlab(da,viewer);
114: #endif
115:   }
116:   return(0);
117: }

119: /*
120:       M is number of grid points
121:       m is number of processors

123: */
126: PetscErrorCode  DMDASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm)
127: {
129:   PetscInt       m,n = 0,x = 0,y = 0;
130:   PetscMPIInt    size,csize,rank;

133:   MPI_Comm_size(comm,&size);
134:   MPI_Comm_rank(comm,&rank);

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;

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;}

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;

157:     MPI_Comm_group(comm,&entire_group);
158:     PetscMalloc1(csize,&groupies);
159:     for (i=0; i<csize; i++) {
160:       groupies[i] = (rank/csize)*csize + i;
161:     }
162:     MPI_Group_incl(entire_group,csize,groupies,&sub_group);
163:     PetscFree(groupies);
164:     MPI_Comm_create(comm,sub_group,outcomm);
165:     MPI_Group_free(&entire_group);
166:     MPI_Group_free(&sub_group);
167:     PetscInfo1(0,"DMDASplitComm2d:Creating redundant coarse problems of size %d\n",csize);
168:   } else {
169:     *outcomm = comm;
170:   }
171:   return(0);
172: }

174: #if defined(new)
177: /*
178:   DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
179:     function lives on a DMDA

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;
191:   PetscInt       gI,nI;
192:   MatStencil     stencil;
193:   DMDALocalInfo  info;

196:   (*ctx->func)(0,U,a,ctx->funcctx);
197:   (*ctx->funcisetbase)(U,ctx->funcctx);

199:   VecGetArray(U,&ww);
200:   VecGetArray(a,&aa);

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;

209:   ww[gI] += h;
210:   (*ctx->funci)(i,w,&v,ctx->funcctx);
211:   aa[nI]  = (v - aa[nI])/h;
212:   ww[gI] -= h;
213:   nI++;

215:   VecRestoreArray(U,&ww);
216:   VecRestoreArray(a,&aa);
217:   return(0);
218: }
219: #endif

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:   DMDABoundaryType bx           = dd->bx;
233:   DMDABoundaryType 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,ltogis;
247:   PetscErrorCode   ierr;

250:   if (stencil_type == DMDA_STENCIL_BOX && (bx == DMDA_BOUNDARY_MIRROR || by == DMDA_BOUNDARY_MIRROR)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Mirror boundary and box stencil");
251:   PetscObjectGetComm((PetscObject)da,&comm);
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

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);

259:   MPI_Comm_size(comm,&size);
260:   MPI_Comm_rank(comm,&rank);

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:   }

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");

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);

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:     PetscMalloc1(m, &dd->lx);
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

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:     PetscMalloc1(n, &dd->ly);
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

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 == 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);
346:   if ((y < s) && ((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);
347:   xe = xs + x;
348:   ye = ys + y;

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:   }

372:   if (bx == DMDA_BOUNDARY_PERIODIC || bx == DMDA_BOUNDARY_MIRROR) {
373:     IXs = xs - s;
374:     IXe = xe + s;
375:     Xs  = xs - s;
376:     Xe  = xe + s;
377:   }

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:   }

400:   if (by == DMDA_BOUNDARY_PERIODIC || by == DMDA_BOUNDARY_MIRROR) {
401:     IYs = ys - s;
402:     IYe = ye + s;
403:     Ys  = ys - s;
404:     Ye  = ye + s;
405:   }

407:   /* stencil length in each direction */
408:   s_x = s;
409:   s_y = s;

411:   /* determine starting point of each processor */
412:   nn       = x*y;
413:   PetscMalloc2(size+1,&bases,size,&ldims);
414:   MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);
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;

424:   /* allocate the base parallel and sequential vectors */
425:   dd->Nlocal = x*y*dof;
426:   VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,0,&global);
427:   dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof;
428:   VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,0,&local);

430:   /* generate appropriate vector scatters */
431:   /* local to global inserts non-ghost point region into global */
432:   VecGetOwnershipRange(global,&start,&end);
433:   ISCreateStride(comm,x*y*dof,start,1,&to);

435:   PetscMalloc1(x*y,&idx);
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:   }

445:   ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);
446:   VecScatterCreate(local,from,global,to,&ltog);
447:   PetscLogObjectParent((PetscObject)da,(PetscObject)ltog);
448:   ISDestroy(&from);
449:   ISDestroy(&to);

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:     count = (IXe-IXs)*(IYe-IYs);
455:     PetscMalloc1(count,&idx);

457:     left  = IXs - Xs; right = left + (IXe-IXs);
458:     down  = IYs - Ys; up = down + (IYe-IYs);
459:     count = 0;
460:     for (i=down; i<up; i++) {
461:       for (j=left; j<right; j++) {
462:         idx[count++] = j + i*(Xe-Xs);
463:       }
464:     }
465:     ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);

467:   } else {
468:     /* must drop into cross shape region */
469:     /*       ---------|
470:             |  top    |
471:          |---         ---| up
472:          |   middle      |
473:          |               |
474:          ----         ---- down
475:             | bottom  |
476:             -----------
477:          Xs xs        xe Xe */
478:     count = (ys-IYs)*x + y*(IXe-IXs) + (IYe-ye)*x;
479:     PetscMalloc1(count,&idx);

481:     left  = xs - Xs; right = left + x;
482:     down  = ys - Ys; up = down + y;
483:     count = 0;
484:     /* bottom */
485:     for (i=(IYs-Ys); i<down; i++) {
486:       for (j=left; j<right; j++) {
487:         idx[count++] = j + i*(Xe-Xs);
488:       }
489:     }
490:     /* middle */
491:     for (i=down; i<up; i++) {
492:       for (j=(IXs-Xs); j<(IXe-Xs); j++) {
493:         idx[count++] = j + i*(Xe-Xs);
494:       }
495:     }
496:     /* top */
497:     for (i=up; i<up+IYe-ye; i++) {
498:       for (j=left; j<right; j++) {
499:         idx[count++] = j + i*(Xe-Xs);
500:       }
501:     }
502:     ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);
503:   }


506:   /* determine who lies on each side of us stored in    n6 n7 n8
507:                                                         n3    n5
508:                                                         n0 n1 n2
509:   */

511:   /* Assume the Non-Periodic Case */
512:   n1 = rank - m;
513:   if (rank % m) {
514:     n0 = n1 - 1;
515:   } else {
516:     n0 = -1;
517:   }
518:   if ((rank+1) % m) {
519:     n2 = n1 + 1;
520:     n5 = rank + 1;
521:     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
522:   } else {
523:     n2 = -1; n5 = -1; n8 = -1;
524:   }
525:   if (rank % m) {
526:     n3 = rank - 1;
527:     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
528:   } else {
529:     n3 = -1; n6 = -1;
530:   }
531:   n7 = rank + m; if (n7 >= m*n) n7 = -1;

533:   if (bx == DMDA_BOUNDARY_PERIODIC && by == DMDA_BOUNDARY_PERIODIC) {
534:     /* Modify for Periodic Cases */
535:     /* Handle all four corners */
536:     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
537:     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
538:     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
539:     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;

541:     /* Handle Top and Bottom Sides */
542:     if (n1 < 0) n1 = rank + m * (n-1);
543:     if (n7 < 0) n7 = rank - m * (n-1);
544:     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
545:     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
546:     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
547:     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;

549:     /* Handle Left and Right Sides */
550:     if (n3 < 0) n3 = rank + (m-1);
551:     if (n5 < 0) n5 = rank - (m-1);
552:     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
553:     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
554:     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
555:     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
556:   } else if (by == DMDA_BOUNDARY_PERIODIC) {  /* Handle Top and Bottom Sides */
557:     if (n1 < 0) n1 = rank + m * (n-1);
558:     if (n7 < 0) n7 = rank - m * (n-1);
559:     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
560:     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
561:     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
562:     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
563:   } else if (bx == DMDA_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */
564:     if (n3 < 0) n3 = rank + (m-1);
565:     if (n5 < 0) n5 = rank - (m-1);
566:     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
567:     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
568:     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
569:     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
570:   }

572:   PetscMalloc1(9,&dd->neighbors);

574:   dd->neighbors[0] = n0;
575:   dd->neighbors[1] = n1;
576:   dd->neighbors[2] = n2;
577:   dd->neighbors[3] = n3;
578:   dd->neighbors[4] = rank;
579:   dd->neighbors[5] = n5;
580:   dd->neighbors[6] = n6;
581:   dd->neighbors[7] = n7;
582:   dd->neighbors[8] = n8;

584:   if (stencil_type == DMDA_STENCIL_STAR) {
585:     /* save corner processor numbers */
586:     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
587:     n0  = n2 = n6 = n8 = -1;
588:   }

590:   PetscMalloc1((Xe-Xs)*(Ye-Ys),&idx);
591:   PetscLogObjectMemory((PetscObject)da,(Xe-Xs)*(Ye-Ys)*sizeof(PetscInt));

593:   nn = 0;
594:   xbase = bases[rank];
595:   for (i=1; i<=s_y; i++) {
596:     if (n0 >= 0) { /* left below */
597:       x_t = lx[n0 % m];
598:       y_t = ly[(n0/m)];
599:       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
600:       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
601:     }

603:     if (n1 >= 0) { /* directly below */
604:       x_t = x;
605:       y_t = ly[(n1/m)];
606:       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
607:       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
608:     } else if (by == DMDA_BOUNDARY_MIRROR) {
609:       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
610:     }

612:     if (n2 >= 0) { /* right below */
613:       x_t = lx[n2 % m];
614:       y_t = ly[(n2/m)];
615:       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
616:       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
617:     }
618:   }

620:   for (i=0; i<y; i++) {
621:     if (n3 >= 0) { /* directly left */
622:       x_t = lx[n3 % m];
623:       /* y_t = y; */
624:       s_t = bases[n3] + (i+1)*x_t - s_x;
625:       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
626:     } else if (bx == DMDA_BOUNDARY_MIRROR) {
627:       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
628:     }

630:     for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */

632:     if (n5 >= 0) { /* directly right */
633:       x_t = lx[n5 % m];
634:       /* y_t = y; */
635:       s_t = bases[n5] + (i)*x_t;
636:       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
637:     } else if (bx == DMDA_BOUNDARY_MIRROR) {
638:       for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
639:     }
640:   }

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:     }

650:     if (n7 >= 0) { /* directly above */
651:       x_t = x;
652:       /* y_t = ly[(n7/m)]; */
653:       s_t = bases[n7] + (i-1)*x_t;
654:       for (j=0; j<x_t; j++) idx[nn++] = s_t++;
655:     } else if (by == DMDA_BOUNDARY_MIRROR) {
656:       for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
657:     }

659:     if (n8 >= 0) { /* right above */
660:       x_t = lx[n8 % m];
661:       /* y_t = ly[(n8/m)]; */
662:       s_t = bases[n8] + (i-1)*x_t;
663:       for (j=0; j<s_x; j++) idx[nn++] = s_t++;
664:     }
665:   }

667:   ISCreateBlock(comm,dof,nn,idx,PETSC_COPY_VALUES,&from);
668:   VecScatterCreate(global,from,local,to,&gtol);
669:   PetscLogObjectParent((PetscObject)da,(PetscObject)gtol);
670:   ISDestroy(&to);
671:   ISDestroy(&from);

673:   if (stencil_type == DMDA_STENCIL_STAR) {
674:     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
675:   }

677:   if (((stencil_type == DMDA_STENCIL_STAR)  ||
678:        (bx && bx != DMDA_BOUNDARY_PERIODIC) ||
679:        (by && by != DMDA_BOUNDARY_PERIODIC))) {
680:     /*
681:         Recompute the local to global mappings, this time keeping the
682:       information about the cross corner processor numbers and any ghosted
683:       but not periodic indices.
684:     */
685:     nn    = 0;
686:     xbase = bases[rank];
687:     for (i=1; i<=s_y; i++) {
688:       if (n0 >= 0) { /* left below */
689:         x_t = lx[n0 % m];
690:         y_t = ly[(n0/m)];
691:         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
692:         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
693:       } else if (xs-Xs > 0 && ys-Ys > 0) {
694:         for (j=0; j<s_x; j++) idx[nn++] = -1;
695:       }
696:       if (n1 >= 0) { /* directly below */
697:         x_t = x;
698:         y_t = ly[(n1/m)];
699:         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
700:         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
701:       } else if (ys-Ys > 0) {
702:         if (by == DMDA_BOUNDARY_MIRROR) {
703:           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(s_y - i + 1)  + j;
704:         } else {
705:           for (j=0; j<x; j++) idx[nn++] = -1;
706:         }
707:       }
708:       if (n2 >= 0) { /* right below */
709:         x_t = lx[n2 % m];
710:         y_t = ly[(n2/m)];
711:         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
712:         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
713:       } else if (Xe-xe> 0 && ys-Ys > 0) {
714:         for (j=0; j<s_x; j++) idx[nn++] = -1;
715:       }
716:     }

718:     for (i=0; i<y; i++) {
719:       if (n3 >= 0) { /* directly left */
720:         x_t = lx[n3 % m];
721:         /* y_t = y; */
722:         s_t = bases[n3] + (i+1)*x_t - s_x;
723:         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
724:       } else if (xs-Xs > 0) {
725:         if (bx == DMDA_BOUNDARY_MIRROR) {
726:           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*i + s_x - j;
727:         } else {
728:           for (j=0; j<s_x; j++) idx[nn++] = -1;
729:         }
730:       }

732:       for (j=0; j<x; j++) idx[nn++] = xbase++; /* interior */

734:       if (n5 >= 0) { /* directly right */
735:         x_t = lx[n5 % m];
736:         /* y_t = y; */
737:         s_t = bases[n5] + (i)*x_t;
738:         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
739:       } else if (Xe-xe > 0) {
740:         if (bx == DMDA_BOUNDARY_MIRROR) {
741:           for (j=0; j<s_x; j++) idx[nn++] = bases[rank] + x*(i + 1) - 2 - j;
742:         } else {
743:           for (j=0; j<s_x; j++) idx[nn++] = -1;
744:         }
745:       }
746:     }

748:     for (i=1; i<=s_y; i++) {
749:       if (n6 >= 0) { /* left above */
750:         x_t = lx[n6 % m];
751:         /* y_t = ly[(n6/m)]; */
752:         s_t = bases[n6] + (i)*x_t - s_x;
753:         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
754:       } else if (xs-Xs > 0 && Ye-ye > 0) {
755:         for (j=0; j<s_x; j++) idx[nn++] = -1;
756:       }
757:       if (n7 >= 0) { /* directly above */
758:         x_t = x;
759:         /* y_t = ly[(n7/m)]; */
760:         s_t = bases[n7] + (i-1)*x_t;
761:         for (j=0; j<x_t; j++) idx[nn++] = s_t++;
762:       } else if (Ye-ye > 0) {
763:         if (by == DMDA_BOUNDARY_MIRROR) {
764:           for (j=0; j<x; j++) idx[nn++] = bases[rank] + x*(y - i - 1)  + j;
765:         } else {
766:           for (j=0; j<x; j++) idx[nn++] = -1;
767:         }
768:       }
769:       if (n8 >= 0) { /* right above */
770:         x_t = lx[n8 % m];
771:         /* y_t = ly[(n8/m)]; */
772:         s_t = bases[n8] + (i-1)*x_t;
773:         for (j=0; j<s_x; j++) idx[nn++] = s_t++;
774:       } else if (Xe-xe > 0 && Ye-ye > 0) {
775:         for (j=0; j<s_x; j++) idx[nn++] = -1;
776:       }
777:     }
778:   }
779:   /*
780:      Set the local to global ordering in the global vector, this allows use
781:      of VecSetValuesLocal().
782:   */
783:   ISCreateBlock(comm,dof,nn,idx,PETSC_OWN_POINTER,&ltogis);
784:   ISLocalToGlobalMappingCreateIS(ltogis,&da->ltogmap);
785:   PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);
786:   ISDestroy(&ltogis);
787:   ISLocalToGlobalMappingBlock(da->ltogmap,dd->w,&da->ltogmapb);
788:   PetscLogObjectParent((PetscObject)da,(PetscObject)da->ltogmap);

790:   PetscFree2(bases,ldims);
791:   dd->m = m;  dd->n  = n;
792:   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
793:   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
794:   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;

796:   VecDestroy(&local);
797:   VecDestroy(&global);

799:   dd->gtol      = gtol;
800:   dd->ltog      = ltog;
801:   dd->base      = base;
802:   da->ops->view = DMView_DA_2d;
803:   dd->ltol      = NULL;
804:   dd->ao        = NULL;
805:   return(0);
806: }

810: /*@C
811:    DMDACreate2d -  Creates an object that will manage the communication of  two-dimensional
812:    regular array data that is distributed across some processors.

814:    Collective on MPI_Comm

816:    Input Parameters:
817: +  comm - MPI communicator
818: .  bx,by - type of ghost nodes the array have.
819:          Use one of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC.
820: .  stencil_type - stencil type.  Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
821: .  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
822:             from the command line with -da_grid_x <M> -da_grid_y <N>)
823: .  m,n - corresponding number of processors in each dimension
824:          (or PETSC_DECIDE to have calculated)
825: .  dof - number of degrees of freedom per node
826: .  s - stencil width
827: -  lx, ly - arrays containing the number of nodes in each cell along
828:            the x and y coordinates, or NULL. If non-null, these
829:            must be of length as m and n, and the corresponding
830:            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
831:            must be M, and the sum of the ly[] entries must be N.

833:    Output Parameter:
834: .  da - the resulting distributed array object

836:    Options Database Key:
837: +  -dm_view - Calls DMView() at the conclusion of DMDACreate2d()
838: .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
839: .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
840: .  -da_processors_x <nx> - number of processors in x direction
841: .  -da_processors_y <ny> - number of processors in y direction
842: .  -da_refine_x <rx> - refinement ratio in x direction
843: .  -da_refine_y <ry> - refinement ratio in y direction
844: -  -da_refine <n> - refine the DMDA n times before creating, if M or N < 0


847:    Level: beginner

849:    Notes:
850:    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
851:    standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
852:    the standard 9-pt stencil.

854:    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
855:    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
856:    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.

858: .keywords: distributed array, create, two-dimensional

860: .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
861:           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(),
862:           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()

864: @*/

866: PetscErrorCode  DMDACreate2d(MPI_Comm comm,DMDABoundaryType bx,DMDABoundaryType by,DMDAStencilType stencil_type,
867:                           PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
868: {

872:   DMDACreate(comm, da);
873:   DMDASetDim(*da, 2);
874:   DMDASetSizes(*da, M, N, 1);
875:   DMDASetNumProcs(*da, m, n, PETSC_DECIDE);
876:   DMDASetBoundaryType(*da, bx, by, DMDA_BOUNDARY_NONE);
877:   DMDASetDof(*da, dof);
878:   DMDASetStencilType(*da, stencil_type);
879:   DMDASetStencilWidth(*da, s);
880:   DMDASetOwnershipRanges(*da, lx, ly, NULL);
881:   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
882:   DMSetFromOptions(*da);
883:   DMSetUp(*da);
884:   return(0);
885: }