Actual source code: da2.c

  1: /*$Id: da2.c,v 1.164 2001/04/10 19:37:23 bsmith Exp bsmith $*/
  2: 
 3:  #include src/dm/da/daimpl.h

  5: int DAGetOwnershipRange(DA da,int **lx,int **ly,int **lz)
  6: {
  9:   if (lx) *lx = da->lx;
 10:   if (ly) *ly = da->ly;
 11:   if (lz) *lz = da->lz;
 12:   return(0);
 13: }

 15: int DAView_2d(DA da,PetscViewer viewer)
 16: {
 17:   int        rank,ierr;
 18:   PetscTruth isascii,isdraw,isbinary;

 21:   MPI_Comm_rank(da->comm,&rank);

 23:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
 24:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
 25:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
 26:   if (isascii) {
 27:     PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %d N %d m %d n %d w %d s %dn",rank,da->M,
 28:                              da->N,da->m,da->n,da->w,da->s);
 29:     PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %d %d, Y range of indices: %d %dn",da->xs,da->xe,da->ys,da->ye);
 30:     PetscViewerFlush(viewer);
 31:   } else if (isdraw) {
 32:     PetscDraw       draw;
 33:     double     ymin = -1*da->s-1,ymax = da->N+da->s;
 34:     double     xmin = -1*da->s-1,xmax = da->M+da->s;
 35:     double     x,y;
 36:     int        base,*idx;
 37:     char       node[10];
 38:     PetscTruth isnull;
 39: 
 40:     PetscViewerDrawGetDraw(viewer,0,&draw);
 41:     PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
 42:     PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);
 43:     PetscDrawSynchronizedClear(draw);

 45:     /* first processor draw all node lines */
 46:     if (!rank) {
 47:       ymin = 0.0; ymax = da->N - 1;
 48:       for (xmin=0; xmin<da->M; xmin++) {
 49:         PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);
 50:       }
 51:       xmin = 0.0; xmax = da->M - 1;
 52:       for (ymin=0; ymin<da->N; ymin++) {
 53:         PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);
 54:       }
 55:     }
 56:     PetscDrawSynchronizedFlush(draw);
 57:     PetscDrawPause(draw);

 59:     /* draw my box */
 60:     ymin = da->ys; ymax = da->ye - 1; xmin = da->xs/da->w;
 61:     xmax =(da->xe-1)/da->w;
 62:     PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);
 63:     PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);
 64:     PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);
 65:     PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);

 67:     /* put in numbers */
 68:     base = (da->base)/da->w;
 69:     for (y=ymin; y<=ymax; y++) {
 70:       for (x=xmin; x<=xmax; x++) {
 71:         sprintf(node,"%d",base++);
 72:         PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);
 73:       }
 74:     }

 76:     PetscDrawSynchronizedFlush(draw);
 77:     PetscDrawPause(draw);
 78:     /* overlay ghost numbers, useful for error checking */
 79:     /* put in numbers */

 81:     base = 0; idx = da->idx;
 82:     ymin = da->Ys; ymax = da->Ye; xmin = da->Xs; xmax = da->Xe;
 83:     for (y=ymin; y<ymax; y++) {
 84:       for (x=xmin; x<xmax; x++) {
 85:         if ((base % da->w) == 0) {
 86:           sprintf(node,"%d",idx[base]/da->w);
 87:           PetscDrawString(draw,x/da->w,y,PETSC_DRAW_BLUE,node);
 88:         }
 89:         base++;
 90:       }
 91:     }
 92:     PetscDrawSynchronizedFlush(draw);
 93:     PetscDrawPause(draw);
 94:   } else if (isbinary) {
 95:     DAView_Binary(da,viewer);
 96:   } else {
 97:     SETERRQ1(1,"Viewer type %s not supported for DA2d",((PetscObject)viewer)->type_name);
 98:   }
 99:   return(0);
100: }

102: #if defined(PETSC_HAVE_AMS)
103: /*
104:       This function tells the AMS the layout of the vectors, it is called
105:    in the VecPublish_xx routines.
106: */
107: EXTERN_C_BEGIN
108: int AMSSetFieldBlock_DA(AMS_Memory amem,char *name,Vec vec)
109: {
110:   int        ierr,dof,dim,ends[4],shift = 0,starts[] = {0,0,0,0};
111:   DA         da = 0;
112:   PetscTruth isseq,ismpi;

115:   if (((PetscObject)vec)->amem < 0) return(0); /* return if not published */

117:   PetscObjectQuery((PetscObject)vec,"DA",(PetscObject*)&da);
118:   if (!da) return(0);
119:   DAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0);
120:   if (dof > 1) {dim++; shift = 1; ends[0] = dof;}

122:   PetscTypeCompare((PetscObject)vec,VEC_SEQ,&isseq);
123:   PetscTypeCompare((PetscObject)vec,VEC_MPI,&ismpi);
124:   if (isseq) {
125:     DAGetGhostCorners(da,0,0,0,ends+shift,ends+shift+1,ends+shift+2);
126:     ends[shift]   += starts[shift]-1;
127:     ends[shift+1] += starts[shift+1]-1;
128:     ends[shift+2] += starts[shift+2]-1;
129:     AMS_Memory_set_field_block(amem,name,dim,starts,ends);
130:     if (ierr) {
131:       char *message;
132:       AMS_Explain_error(ierr,&message);
133:       SETERRQ(ierr,message);
134:     }
135:   } else if (ismpi) {
136:     DAGetCorners(da,starts+shift,starts+shift+1,starts+shift+2,
137:                            ends+shift,ends+shift+1,ends+shift+2);
138:     ends[shift]   += starts[shift]-1;
139:     ends[shift+1] += starts[shift+1]-1;
140:     ends[shift+2] += starts[shift+2]-1;
141:     AMS_Memory_set_field_block(amem,name,dim,starts,ends);
142:     if (ierr) {
143:       char *message;
144:       AMS_Explain_error(ierr,&message);
145:       SETERRQ(ierr,message);
146:     }
147:   } else {
148:     SETERRQ1(1,"Wrong vector type %s for this call",((PetscObject)vec)->type_name);
149:   }

151:   return(0);
152: }
153: EXTERN_C_END
154: #endif

156: int DAPublish_Petsc(PetscObject obj)
157: {
158: #if defined(PETSC_HAVE_AMS)
159:   DA          v = (DA) obj;
160:   int         ierr;
161: #endif


165: #if defined(PETSC_HAVE_AMS)
166:   /* if it is already published then return */
167:   if (v->amem >=0) return(0);

169:   PetscObjectPublishBaseBegin(obj);
170:   PetscObjectPublishBaseEnd(obj);
171: #endif

173:   return(0);
174: }

176: /*
177:    This allows the DA vectors to properly tell Matlab their dimensions
178: */
179: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX)
180: #include "engine.h"   /* Matlab include file */
181: #include "mex.h"      /* Matlab include file */
182: EXTERN_C_BEGIN
183: int VecMatlabEnginePut_DA2d(PetscObject obj,void *engine)
184: {
185:   int     ierr,n,m;
186:   Vec     vec = (Vec)obj;
187:   Scalar  *array;
188:   mxArray *mat;
189:   DA      da;

192:   PetscObjectQuery((PetscObject)vec,"DA",(PetscObject*)&da);
193:   if (!da) SETERRQ(1,"Vector not associated with a DA");
194:   DAGetGhostCorners(da,0,0,0,&m,&n,0);

196:   VecGetArray(vec,&array);
197: #if !defined(PETSC_USE_COMPLEX)
198:   mat  = mxCreateDoubleMatrix(m,n,mxREAL);
199: #else
200:   mat  = mxCreateDoubleMatrix(m,n,mxCOMPLEX);
201: #endif
202:   PetscMemcpy(mxGetPr(mat),array,n*m*sizeof(Scalar));
203:   PetscObjectName(obj);
204:   mxSetName(mat,obj->name);
205:   engPutArray((Engine *)engine,mat);
206: 
207:   VecRestoreArray(vec,&array);
208:   return(0);
209: }
210: EXTERN_C_END
211: #endif

213: /*@C
214:    DACreate2d -  Creates an object that will manage the communication of  two-dimensional 
215:    regular array data that is distributed across some processors.

217:    Collective on MPI_Comm

219:    Input Parameters:
220: +  comm - MPI communicator
221: .  wrap - type of periodicity should the array have. 
222:          Use one of DA_NONPERIODIC, DA_XPERIODIC, DA_YPERIODIC, or DA_XYPERIODIC.
223: .  stencil_type - stencil type.  Use either DA_STENCIL_BOX or DA_STENCIL_STAR.
224: .  M,N - global dimension in each direction of the array
225: .  m,n - corresponding number of processors in each dimension 
226:          (or PETSC_DECIDE to have calculated)
227: .  dof - number of degrees of freedom per node
228: .  s - stencil width
229: -  lx, ly - arrays containing the number of nodes in each cell along
230:            the x and y coordinates, or PETSC_NULL. If non-null, these
231:            must be of length as m and n, and the corresponding
232:            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
233:            must be M, and the sum of the ly[] entries must be N.

235:    Output Parameter:
236: .  inra - the resulting distributed array object

238:    Options Database Key:
239: +  -da_view - Calls DAView() at the conclusion of DACreate2d()
240: .  -da_grid_x <nx> - number of grid points in x direction
241: .  -da_grid_y <ny> - number of grid points in y direction
242: .  -da_processors_x <nx> - number of processors in x direction
243: .  -da_processors_y <ny> - number of processors in y direction
244: -  -da_noao - do not compute natural to PETSc ordering object

246:    Level: beginner

248:    Notes:
249:    The stencil type DA_STENCIL_STAR with width 1 corresponds to the 
250:    standard 5-pt stencil, while DA_STENCIL_BOX with width 1 denotes
251:    the standard 9-pt stencil.

253:    The array data itself is NOT stored in the DA, it is stored in Vec objects;
254:    The appropriate vector objects can be obtained with calls to DACreateGlobalVector()
255:    and DACreateLocalVector() and calls to VecDuplicate() if more are needed.

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

259: .seealso: DADestroy(), DAView(), DACreate1d(), DACreate3d(), DAGlobalToLocalBegin(),
260:           DAGlobalToLocalEnd(), DALocalToGlobal(), DALocalToLocalBegin(), DALocalToLocalEnd(),
261:           DAGetInfo(), DACreateGlobalVector(), DACreateLocalVector(), DACreateNaturalVector(), DALoad(), DAView()

263: @*/
264: int DACreate2d(MPI_Comm comm,DAPeriodicType wrap,DAStencilType stencil_type,
265:                 int M,int N,int m,int n,int dof,int s,int *lx,int *ly,DA *inra)
266: {
267:   int           rank,size,xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,ierr,start,end;
268:   int           up,down,left,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn;
269:   int           xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
270:   int           s_x,s_y; /* s proportionalized to w */
271:   int           *gA,*gB,*gAall,*gBall,ict,ldim,gdim,*flx = 0,*fly = 0;
272:   int           sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0,refine_x = 2, refine_y = 2;
273:   PetscTruth    flg1,flg2;
274:   DA            da;
275:   Vec           local,global;
276:   VecScatter    ltog,gtol;
277:   IS            to,from;

280:   *inra = 0;

282:   if (dof < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %d",dof);
283:   if (s < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %d",s);
284:   if (M < 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Must have M positive");
285:   if (N < 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Must have N positive");

287:   PetscOptionsBegin(comm,PETSC_NULL,"2d DA Options","DA");
288:     PetscOptionsInt("-da_grid_x","Number of grid points in x direction","DACreate2d",M,&M,PETSC_NULL);
289:     PetscOptionsInt("-da_grid_y","Number of grid points in y direction","DACreate2d",N,&N,PETSC_NULL);
290:     PetscOptionsInt("-da_processors_x","Number of processors in x direction","DACreate2d",m,&m,PETSC_NULL);
291:     PetscOptionsInt("-da_processors_y","Number of processors in y direction","DACreate2d",n,&n,PETSC_NULL);
292:     PetscOptionsInt("-da_refine_x","Refinement ratio in x direction","DACreate2d",refine_x,&refine_x,PETSC_NULL);
293:     PetscOptionsInt("-da_refine_y","Refinement ratio in y direction","DACreate2d",refine_y,&refine_y,PETSC_NULL);
294:   PetscOptionsEnd();

296:   PetscHeaderCreate(da,_p_DA,struct _DAOps,DA_COOKIE,0,"DA",comm,DADestroy,DAView);
297:   PetscLogObjectCreate(da);
298:   da->bops->publish           = DAPublish_Petsc;
299:   da->ops->createglobalvector = DACreateGlobalVector;
300:   da->ops->getinterpolation   = DAGetInterpolation;
301:   da->ops->getcoloring        = DAGetColoring;
302:   da->ops->refine             = DARefine;
303:   PetscLogObjectMemory(da,sizeof(struct _p_DA));
304:   da->dim        = 2;
305:   da->gtog1      = 0;
306:   da->refine_x   = refine_x;
307:   da->refine_y   = refine_y;
308:   PetscMalloc(dof*sizeof(char*),&da->fieldname);
309:   PetscMemzero(da->fieldname,dof*sizeof(char*));

311:   MPI_Comm_size(comm,&size);
312:   MPI_Comm_rank(comm,&rank);

314:   if (m != PETSC_DECIDE) {
315:     if (m < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %d",m);}
316:     else if (m > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %d %d",m,size);}
317:   }
318:   if (n != PETSC_DECIDE) {
319:     if (n < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %d",n);}
320:     else if (n > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %d %d",n,size);}
321:   }

323:   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
324:     /* try for squarish distribution */
325:     /* This should use MPI_Dims_create instead */
326:     m = (int)(0.5 + sqrt(((double)M)*((double)size)/((double)N)));
327:     if (!m) m = 1;
328:     while (m > 0) {
329:       n = size/m;
330:       if (m*n == size) break;
331:       m--;
332:     }
333:     if (M > N && m < n) {int _m = m; m = n; n = _m;}
334:     if (m*n != size) SETERRQ(PETSC_ERR_PLIB,"Internally Created Bad Partition");
335:   } else if (m*n != size) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");

337:   if (M < m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %d %d",M,m);
338:   if (N < n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %d %d",N,n);

340:   /*
341:      We should create an MPI Cartesian topology here, with reorder
342:      set to true.  That would create a NEW communicator that we would
343:      need to use for operations on this distributed array 
344:   */
345:   PetscOptionsHasName(PETSC_NULL,"-da_partition_nodes_at_end",&flg2);

347:   /* 
348:      Determine locally owned region 
349:      xs is the first local node number, x is the number of local nodes 
350:   */
351:   if (lx) { /* user sets distribution */
352:     x  = lx[rank % m];
353:     xs = 0;
354:     for (i=0; i<(rank % m); i++) {
355:       xs += lx[i];
356:     }
357:     left = xs;
358:     for (i=(rank % m); i<m; i++) {
359:       left += lx[i];
360:     }
361:     if (left != M) {
362:       SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M: %d %d",left,M);
363:     }
364:   } else if (flg2) {
365:     x = (M + rank%m)/m;
366:     if (x < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column width is too thin for stencil! %d %d",x,s);
367:     if (M/m == x) { xs = (rank % m)*x; }
368:     else          { xs = (rank % m)*(x-1) + (M+(rank % m))%(x*m); }
369:     SETERRQ(PETSC_ERR_SUP,"-da_partition_nodes_at_end not supported");
370:   } else { /* Normal PETSc distribution */
371:     x = M/m + ((M % m) > (rank % m));
372:     if (x < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column width is too thin for stencil! %d %d",x,s);
373:     if ((M % m) > (rank % m)) { xs = (rank % m)*x; }
374:     else                      { xs = (M % m)*(x+1) + ((rank % m)-(M % m))*x; }
375:     PetscMalloc(m*sizeof(int),&lx);
376:     flx = lx;
377:     for (i=0; i<m; i++) {
378:       lx[i] = M/m + ((M % m) > i);
379:     }
380:   }

382:   /* 
383:      Determine locally owned region 
384:      ys is the first local node number, y is the number of local nodes 
385:   */
386:   if (ly) { /* user sets distribution */
387:     y  = ly[rank/m];
388:     ys = 0;
389:     for (i=0; i<(rank/m); i++) {
390:       ys += ly[i];
391:     }
392:     left = ys;
393:     for (i=(rank/m); i<n; i++) {
394:       left += ly[i];
395:     }
396:     if (left != N) {
397:       SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Sum of ly across processors not equal to N: %d %d",left,N);
398:     }
399:   } else if (flg2) {
400:     y = (N + rank/m)/n;
401:     if (y < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row width is too thin for stencil! %d %d",y,s);
402:     if (N/n == y) { ys = (rank/m)*y;  }
403:     else          { ys = (rank/m)*(y-1) + (N+(rank/m))%(y*n); }
404:     SETERRQ(PETSC_ERR_SUP,"-da_partition_nodes_at_end not supported");
405:   } else { /* Normal PETSc distribution */
406:     y = N/n + ((N % n) > (rank/m));
407:     if (y < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row width is too thin for stencil! %d %d",y,s);
408:     if ((N % n) > (rank/m)) { ys = (rank/m)*y; }
409:     else                    { ys = (N % n)*(y+1) + ((rank/m)-(N % n))*y; }
410:     PetscMalloc(n*sizeof(int),&ly);
411:     fly  = ly;
412:     for (i=0; i<n; i++) {
413:       ly[i] = N/n + ((N % n) > i);
414:     }
415:   }

417:   xe = xs + x;
418:   ye = ys + y;

420:   /* determine ghost region */
421:   /* Assume No Periodicity */
422:   if (xs-s > 0) Xs = xs - s; else Xs = 0;
423:   if (ys-s > 0) Ys = ys - s; else Ys = 0;
424:   if (xe+s <= M) Xe = xe + s; else Xe = M;
425:   if (ye+s <= N) Ye = ye + s; else Ye = N;

427:   /* X Periodic */
428:   if (DAXPeriodic(wrap)){
429:     Xs = xs - s;
430:     Xe = xe + s;
431:   }

433:   /* Y Periodic */
434:   if (DAYPeriodic(wrap)){
435:     Ys = ys - s;
436:     Ye = ye + s;
437:   }

439:   /* Resize all X parameters to reflect w */
440:   x   *= dof;
441:   xs  *= dof;
442:   xe  *= dof;
443:   Xs  *= dof;
444:   Xe  *= dof;
445:   s_x = s*dof;
446:   s_y = s;

448:   /* determine starting point of each processor */
449:   nn = x*y;
450:   PetscMalloc((2*size+1)*sizeof(int),&bases);
451:   ldims = (int*)(bases+size+1);
452:   MPI_Allgather(&nn,1,MPI_INT,ldims,1,MPI_INT,comm);
453:   bases[0] = 0;
454:   for (i=1; i<=size; i++) {
455:     bases[i] = ldims[i-1];
456:   }
457:   for (i=1; i<=size; i++) {
458:     bases[i] += bases[i-1];
459:   }

461:   /* allocate the base parallel and sequential vectors */
462:   VecCreateMPI(comm,x*y,PETSC_DECIDE,&global);
463:   VecSetBlockSize(global,dof);
464:   VecCreateSeq(PETSC_COMM_SELF,(Xe-Xs)*(Ye-Ys),&local);
465:   VecSetBlockSize(local,dof);


468:   /* generate appropriate vector scatters */
469:   /* local to global inserts non-ghost point region into global */
470:   VecGetOwnershipRange(global,&start,&end);
471:   ISCreateStride(comm,x*y,start,1,&to);

473:   left  = xs - Xs; down  = ys - Ys; up    = down + y;
474:   PetscMalloc(x*(up - down)*sizeof(int),&idx);
475:   count = 0;
476:   for (i=down; i<up; i++) {
477:     for (j=0; j<x; j++) {
478:       idx[count++] = left + i*(Xe-Xs) + j;
479:     }
480:   }
481:   ISCreateGeneral(comm,count,idx,&from);
482:   PetscFree(idx);

484:   VecScatterCreate(local,from,global,to,&ltog);
485:   PetscLogObjectParent(da,to);
486:   PetscLogObjectParent(da,from);
487:   PetscLogObjectParent(da,ltog);
488:   ISDestroy(from);
489:   ISDestroy(to);

491:   /* global to local must include ghost points */
492:   if (stencil_type == DA_STENCIL_BOX) {
493:     ISCreateStride(comm,(Xe-Xs)*(Ye-Ys),0,1,&to);
494:   } else {
495:     /* must drop into cross shape region */
496:     /*       ---------|
497:             |  top    |
498:          |---         ---|
499:          |   middle      |
500:          |               |
501:          ----         ----
502:             | bottom  |
503:             -----------
504:         Xs xs        xe  Xe */
505:     /* bottom */
506:     left  = xs - Xs; down = ys - Ys; up    = down + y;
507:     count = down*(xe-xs) + (up-down)*(Xe-Xs) + (Ye-Ys-up)*(xe-xs);
508:     ierr  = PetscMalloc(count*sizeof(int),&idx);
509:     count = 0;
510:     for (i=0; i<down; i++) {
511:       for (j=0; j<xe-xs; j++) {
512:         idx[count++] = left + i*(Xe-Xs) + j;
513:       }
514:     }
515:     /* middle */
516:     for (i=down; i<up; i++) {
517:       for (j=0; j<Xe-Xs; j++) {
518:         idx[count++] = i*(Xe-Xs) + j;
519:       }
520:     }
521:     /* top */
522:     for (i=up; i<Ye-Ys; i++) {
523:       for (j=0; j<xe-xs; j++) {
524:         idx[count++] = left + i*(Xe-Xs) + j;
525:       }
526:     }
527:     ISCreateGeneral(comm,count,idx,&to);
528:     PetscFree(idx);
529:   }


532:   /* determine who lies on each side of use stored in    n6 n7 n8
533:                                                          n3    n5
534:                                                          n0 n1 n2
535:   */

537:   /* Assume the Non-Periodic Case */
538:   n1 = rank - m;
539:   if (rank % m) {
540:     n0 = n1 - 1;
541:   } else {
542:     n0 = -1;
543:   }
544:   if ((rank+1) % m) {
545:     n2 = n1 + 1;
546:     n5 = rank + 1;
547:     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
548:   } else {
549:     n2 = -1; n5 = -1; n8 = -1;
550:   }
551:   if (rank % m) {
552:     n3 = rank - 1;
553:     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
554:   } else {
555:     n3 = -1; n6 = -1;
556:   }
557:   n7 = rank + m; if (n7 >= m*n) n7 = -1;


560:   /* Modify for Periodic Cases */
561:   if (wrap == DA_YPERIODIC) {  /* Handle Top and Bottom Sides */
562:     if (n1 < 0) n1 = rank + m * (n-1);
563:     if (n7 < 0) n7 = rank - m * (n-1);
564:     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
565:     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
566:     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
567:     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
568:   } else if (wrap == DA_XPERIODIC) { /* Handle Left and Right Sides */
569:     if (n3 < 0) n3 = rank + (m-1);
570:     if (n5 < 0) n5 = rank - (m-1);
571:     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
572:     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
573:     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
574:     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
575:   } else if (wrap == DA_XYPERIODIC) {

577:     /* Handle all four corners */
578:     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
579:     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
580:     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
581:     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;

583:     /* Handle Top and Bottom Sides */
584:     if (n1 < 0) n1 = rank + m * (n-1);
585:     if (n7 < 0) n7 = rank - m * (n-1);
586:     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
587:     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
588:     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
589:     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;

591:     /* Handle Left and Right Sides */
592:     if (n3 < 0) n3 = rank + (m-1);
593:     if (n5 < 0) n5 = rank - (m-1);
594:     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
595:     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
596:     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
597:     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
598:   }

600:   if (stencil_type == DA_STENCIL_STAR) {
601:     /* save corner processor numbers */
602:     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
603:     n0 = n2 = n6 = n8 = -1;
604:   }

606:   PetscMalloc((x+2*s_x)*(y+2*s_y)*sizeof(int),&idx);
607:   PetscLogObjectMemory(da,(x+2*s_x)*(y+2*s_y)*sizeof(int));
608:   nn = 0;

610:   xbase = bases[rank];
611:   for (i=1; i<=s_y; i++) {
612:     if (n0 >= 0) { /* left below */
613:       x_t = lx[n0 % m]*dof;
614:       y_t = ly[(n0/m)];
615:       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
616:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
617:     }
618:     if (n1 >= 0) { /* directly below */
619:       x_t = x;
620:       y_t = ly[(n1/m)];
621:       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
622:       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
623:     }
624:     if (n2 >= 0) { /* right below */
625:       x_t = lx[n2 % m]*dof;
626:       y_t = ly[(n2/m)];
627:       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
628:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
629:     }
630:   }

632:   for (i=0; i<y; i++) {
633:     if (n3 >= 0) { /* directly left */
634:       x_t = lx[n3 % m]*dof;
635:       /* y_t = y; */
636:       s_t = bases[n3] + (i+1)*x_t - s_x;
637:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
638:     }

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

642:     if (n5 >= 0) { /* directly right */
643:       x_t = lx[n5 % m]*dof;
644:       /* y_t = y; */
645:       s_t = bases[n5] + (i)*x_t;
646:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
647:     }
648:   }

650:   for (i=1; i<=s_y; i++) {
651:     if (n6 >= 0) { /* left above */
652:       x_t = lx[n6 % m]*dof;
653:       /* y_t = ly[(n6/m)]; */
654:       s_t = bases[n6] + (i)*x_t - s_x;
655:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
656:     }
657:     if (n7 >= 0) { /* directly above */
658:       x_t = x;
659:       /* y_t = ly[(n7/m)]; */
660:       s_t = bases[n7] + (i-1)*x_t;
661:       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
662:     }
663:     if (n8 >= 0) { /* right above */
664:       x_t = lx[n8 % m]*dof;
665:       /* y_t = ly[(n8/m)]; */
666:       s_t = bases[n8] + (i-1)*x_t;
667:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
668:     }
669:   }

671:   base = bases[rank];
672:   ISCreateGeneral(comm,nn,idx,&from);
673:   VecScatterCreate(global,from,local,to,&gtol);
674:   PetscLogObjectParent(da,to);
675:   PetscLogObjectParent(da,from);
676:   PetscLogObjectParent(da,gtol);
677:   ISDestroy(to);
678:   ISDestroy(from);

680:   if (stencil_type == DA_STENCIL_STAR) {
681:     /*
682:         Recompute the local to global mappings, this time keeping the 
683:       information about the cross corner processor numbers.
684:     */
685:     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
686:     nn = 0;
687:     xbase = bases[rank];
688:     for (i=1; i<=s_y; i++) {
689:       if (n0 >= 0) { /* left below */
690:         x_t = lx[n0 % m]*dof;
691:         y_t = ly[(n0/m)];
692:         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
693:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
694:       }
695:       if (n1 >= 0) { /* directly below */
696:         x_t = x;
697:         y_t = ly[(n1/m)];
698:         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
699:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
700:       }
701:       if (n2 >= 0) { /* right below */
702:         x_t = lx[n2 % m]*dof;
703:         y_t = ly[(n2/m)];
704:         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
705:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
706:       }
707:     }

709:     for (i=0; i<y; i++) {
710:       if (n3 >= 0) { /* directly left */
711:         x_t = lx[n3 % m]*dof;
712:         /* y_t = y; */
713:         s_t = bases[n3] + (i+1)*x_t - s_x;
714:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
715:       }

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

719:       if (n5 >= 0) { /* directly right */
720:         x_t = lx[n5 % m]*dof;
721:         /* y_t = y; */
722:         s_t = bases[n5] + (i)*x_t;
723:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
724:       }
725:     }

727:     for (i=1; i<=s_y; i++) {
728:       if (n6 >= 0) { /* left above */
729:         x_t = lx[n6 % m]*dof;
730:         /* y_t = ly[(n6/m)]; */
731:         s_t = bases[n6] + (i)*x_t - s_x;
732:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
733:       }
734:       if (n7 >= 0) { /* directly above */
735:         x_t = x;
736:         /* y_t = ly[(n7/m)]; */
737:         s_t = bases[n7] + (i-1)*x_t;
738:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
739:       }
740:       if (n8 >= 0) { /* right above */
741:         x_t = lx[n8 % m]*dof;
742:         /* y_t = ly[(n8/m)]; */
743:         s_t = bases[n8] + (i-1)*x_t;
744:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
745:       }
746:     }
747:   }

749:   da->M  = M;  da->N  = N;  da->m  = m;  da->n  = n;  da->w = dof;  da->s = s;
750:   da->xs = xs; da->xe = xe; da->ys = ys; da->ye = ye; da->zs = 0; da->ze = 1;
751:   da->Xs = Xs; da->Xe = Xe; da->Ys = Ys; da->Ye = Ye; da->Zs = 0; da->Ze = 1;
752:   da->P  = 1;  da->p  = 1;

754:   PetscLogObjectParent(da,global);
755:   PetscLogObjectParent(da,local);

757:   da->global       = global;
758:   da->local        = local;
759:   da->gtol         = gtol;
760:   da->ltog         = ltog;
761:   da->idx          = idx;
762:   da->Nl           = nn;
763:   da->base         = base;
764:   da->wrap         = wrap;
765:   da->ops->view    = DAView_2d;
766:   da->stencil_type = stencil_type;

768:   /* 
769:      Set the local to global ordering in the global vector, this allows use
770:      of VecSetValuesLocal().
771:   */
772:   ISLocalToGlobalMappingCreate(comm,nn,idx,&da->ltogmap);
773:   VecSetLocalToGlobalMapping(da->global,da->ltogmap);
774:   ISLocalToGlobalMappingBlock(da->ltogmap,da->w,&da->ltogmapb);
775:   VecSetLocalToGlobalMappingBlock(da->global,da->ltogmapb);
776:   PetscLogObjectParent(da,da->ltogmap);

778:   *inra = da;

780:   /* recalculate the idx including missed ghost points */
781:   /* Assume the Non-Periodic Case */
782:   n1 = rank - m;
783:   if (rank % m) {
784:     n0 = n1 - 1;
785:   } else {
786:     n0 = -1;
787:   }
788:   if ((rank+1) % m) {
789:     n2 = n1 + 1;
790:     n5 = rank + 1;
791:     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
792:   } else {
793:     n2 = -1; n5 = -1; n8 = -1;
794:   }
795:   if (rank % m) {
796:     n3 = rank - 1;
797:     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
798:   } else {
799:     n3 = -1; n6 = -1;
800:   }
801:   n7 = rank + m; if (n7 >= m*n) n7 = -1;


804:   /* Modify for Periodic Cases */
805:   if (wrap == DA_YPERIODIC) {  /* Handle Top and Bottom Sides */
806:     if (n1 < 0) n1 = rank + m * (n-1);
807:     if (n7 < 0) n7 = rank - m * (n-1);
808:     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
809:     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
810:     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
811:     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
812:   } else if (wrap == DA_XPERIODIC) { /* Handle Left and Right Sides */
813:     if (n3 < 0) n3 = rank + (m-1);
814:     if (n5 < 0) n5 = rank - (m-1);
815:     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
816:     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
817:     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
818:     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
819:   } else if (wrap == DA_XYPERIODIC) {

821:     /* Handle all four corners */
822:     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
823:     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
824:     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
825:     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;

827:     /* Handle Top and Bottom Sides */
828:     if (n1 < 0) n1 = rank + m * (n-1);
829:     if (n7 < 0) n7 = rank - m * (n-1);
830:     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
831:     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
832:     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
833:     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;

835:     /* Handle Left and Right Sides */
836:     if (n3 < 0) n3 = rank + (m-1);
837:     if (n5 < 0) n5 = rank - (m-1);
838:     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
839:     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
840:     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
841:     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
842:   }

844:   nn = 0;

846:   xbase = bases[rank];
847:   for (i=1; i<=s_y; i++) {
848:     if (n0 >= 0) { /* left below */
849:       x_t = lx[n0 % m]*dof;
850:       y_t = ly[(n0/m)];
851:       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
852:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
853:     }
854:     if (n1 >= 0) { /* directly below */
855:       x_t = x;
856:       y_t = ly[(n1/m)];
857:       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
858:       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
859:     }
860:     if (n2 >= 0) { /* right below */
861:       x_t = lx[n2 % m]*dof;
862:       y_t = ly[(n2/m)];
863:       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
864:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
865:     }
866:   }

868:   for (i=0; i<y; i++) {
869:     if (n3 >= 0) { /* directly left */
870:       x_t = lx[n3 % m]*dof;
871:       /* y_t = y; */
872:       s_t = bases[n3] + (i+1)*x_t - s_x;
873:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
874:     }

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

878:     if (n5 >= 0) { /* directly right */
879:       x_t = lx[n5 % m]*dof;
880:       /* y_t = y; */
881:       s_t = bases[n5] + (i)*x_t;
882:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
883:     }
884:   }

886:   for (i=1; i<=s_y; i++) {
887:     if (n6 >= 0) { /* left above */
888:       x_t = lx[n6 % m]*dof;
889:       /* y_t = ly[(n6/m)]; */
890:       s_t = bases[n6] + (i)*x_t - s_x;
891:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
892:     }
893:     if (n7 >= 0) { /* directly above */
894:       x_t = x;
895:       /* y_t = ly[(n7/m)]; */
896:       s_t = bases[n7] + (i-1)*x_t;
897:       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
898:     }
899:     if (n8 >= 0) { /* right above */
900:       x_t = lx[n8 % m]*dof;
901:       /* y_t = ly[(n8/m)]; */
902:       s_t = bases[n8] + (i-1)*x_t;
903:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
904:     }
905:   }
906:   /* keep bases for use at end of routine */
907:   /* PetscFree(bases); */

909:   /* construct the local to local scatter context */
910:   /* 
911:       We simply remap the values in the from part of 
912:     global to local to read from an array with the ghost values 
913:     rather then from the plan array.
914:   */
915:   VecScatterCopy(gtol,&da->ltol);
916:   PetscLogObjectParent(da,da->ltol);
917:   left  = xs - Xs; down  = ys - Ys; up    = down + y;
918:   PetscMalloc(x*(up - down)*sizeof(int),&idx);
919:   count = 0;
920:   for (i=down; i<up; i++) {
921:     for (j=0; j<x; j++) {
922:       idx[count++] = left + i*(Xe-Xs) + j;
923:     }
924:   }
925:   VecScatterRemap(da->ltol,idx,PETSC_NULL);
926:   PetscFree(idx);

928:   /* 
929:      Build the natural ordering to PETSc ordering mappings.
930:   */
931:   PetscOptionsHasName(PETSC_NULL,"-da_noao",&flg1);
932:   if (!flg1) {
933:     IS  ispetsc,isnatural;
934:     int *lidx,lict = 0,Nlocal = (da->xe-da->xs)*(da->ye-da->ys);

936:     ISCreateStride(comm,Nlocal,da->base,1,&ispetsc);

938:     PetscMalloc(Nlocal*sizeof(int),&lidx);
939:     for (j=ys; j<ye; j++) {
940:       for (i=xs; i<xe; i++) {
941:         /*  global number in natural ordering */
942:         lidx[lict++] = i + j*M*dof;
943:       }
944:     }
945:     ISCreateGeneral(comm,Nlocal,lidx,&isnatural);
946:     PetscFree(lidx);


949:     AOCreateBasicIS(isnatural,ispetsc,&da->ao);
950:     PetscLogObjectParent(da,da->ao);
951:     ISDestroy(ispetsc);
952:     ISDestroy(isnatural);
953:   } else {
954:     da->ao = PETSC_NULL;
955:   }

957:   if (!flx) {
958:     PetscMalloc(m*sizeof(int),&flx);
959:     PetscMemcpy(flx,lx,m*sizeof(int));
960:   }
961:   if (!fly) {
962:     PetscMalloc(n*sizeof(int),&fly);
963:     PetscMemcpy(fly,ly,n*sizeof(int));
964:   }
965:   da->lx = flx;
966:   da->ly = fly;

968:   /*
969:      Note the following will be removed soon. Since the functionality 
970:     is replaced by the above.
971:   */
972:   /* Construct the mapping from current global ordering to global
973:      ordering that would be used if only 1 processor were employed.
974:      This mapping is intended only for internal use by discrete
975:      function and matrix viewers.

977:      Note: At this point, x has already been adjusted for multiple
978:      degrees of freedom per node.
979:    */
980:   ldim = x*y;
981:   VecGetSize(global,&gdim);
982:   PetscMalloc(gdim*sizeof(int),&da->gtog1);
983:   PetscLogObjectMemory(da,gdim*sizeof(int));
984:   PetscMalloc((2*(gdim+ldim))*sizeof(int),&gA);
985:   gB        = (int *)(gA + ldim);
986:   gAall     = (int *)(gB + ldim);
987:   gBall     = (int *)(gAall + gdim);

989:   /* Compute local parts of global orderings */
990:   ict = 0;
991:   for (j=ys; j<ye; j++) {
992:     for (i=xs; i<xe; i++) {
993:       /* gA = global number for 1 proc; gB = current global number */
994:       gA[ict] = i + j*M*dof;
995:       gB[ict] = start + ict;
996:       ict++;
997:     }
998:   }
999:   /* Broadcast the orderings */
1000:   MPI_Allgatherv(gA,ldim,MPI_INT,gAall,ldims,bases,MPI_INT,comm);
1001:   MPI_Allgatherv(gB,ldim,MPI_INT,gBall,ldims,bases,MPI_INT,comm);
1002:   for (i=0; i<gdim; i++) da->gtog1[gBall[i]] = gAall[i];
1003:   PetscFree(gA);
1004:   PetscFree(bases);

1006:   PetscOptionsHasName(PETSC_NULL,"-da_view",&flg1);
1007:   if (flg1) {DAView(da,PETSC_VIEWER_STDOUT_(da->comm));}
1008:   PetscOptionsHasName(PETSC_NULL,"-da_view_draw",&flg1);
1009:   if (flg1) {DAView(da,PETSC_VIEWER_DRAW_(da->comm));}
1010:   PetscOptionsHasName(PETSC_NULL,"-help",&flg1);
1011:   if (flg1) {DAPrintHelp(da);}

1013:   PetscPublishAll(da);
1014: #if defined(PETSC_HAVE_AMS)
1015:   PetscObjectComposeFunctionDynamic((PetscObject)global,"AMSSetFieldBlock_C",
1016:          "AMSSetFieldBlock_DA",AMSSetFieldBlock_DA);
1017:   PetscObjectComposeFunctionDynamic((PetscObject)local,"AMSSetFieldBlock_C",
1018:          "AMSSetFieldBlock_DA",AMSSetFieldBlock_DA);
1019:   if (((PetscObject)global)->amem > -1) {
1020:     AMSSetFieldBlock_DA(((PetscObject)global)->amem,"values",global);
1021:   }
1022: #endif
1023: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX)
1024:   if (dof == 1) {
1025:     PetscObjectComposeFunctionDynamic((PetscObject)local,"PetscMatlabEnginePut_C","VecMatlabEnginePut_DA2d",VecMatlabEnginePut_DA2d);
1026:   }
1027: #endif
1028:   VecSetOperation(global,VECOP_VIEW,(void(*)())VecView_MPI_DA);
1029:   VecSetOperation(global,VECOP_LOADINTOVECTOR,(void(*)())VecLoadIntoVector_Binary_DA);
1030:   return(0);
1031: }

1033: /*@
1034:    DAPrintHelp - Prints command line options for DA.

1036:    Collective on DA

1038:    Input Parameters:
1039: .  da - the distributed array

1041:    Level: intermediate

1043: .seealso: DACreate1d(), DACreate2d(), DACreate3d()

1045: .keywords: DA, help

1047: @*/
1048: int DAPrintHelp(DA da)
1049: {
1050:   static PetscTruth called = PETSC_FALSE;
1051:   MPI_Comm          comm;
1052:   int               ierr;


1057:   comm = da->comm;
1058:   if (!called) {
1059:     (*PetscHelpPrintf)(comm,"General Distributed Array (DA) options:n");
1060:     (*PetscHelpPrintf)(comm,"  -da_view: print DA distribution to screenn");
1061:     (*PetscHelpPrintf)(comm,"  -da_view_draw: display DA in windown");
1062:     called = PETSC_TRUE;
1063:   }
1064:   return(0);
1065: }

1067: /*@C
1068:    DARefine - Creates a new distributed array that is a refinement of a given
1069:    distributed array.

1071:    Collective on DA

1073:    Input Parameter:
1074: +  da - initial distributed array
1075: -  comm - communicator to contain refined DA, must be either same as the da communicator or include the 
1076:           da communicator and be 2, 4, or 8 times larger. Currently ignored

1078:    Output Parameter:
1079: .  daref - refined distributed array

1081:    Level: advanced

1083:    Note:
1084:    Currently, refinement consists of just doubling the number of grid spaces
1085:    in each dimension of the DA.

1087: .keywords:  distributed array, refine

1089: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy()
1090: @*/
1091: int DARefine(DA da,MPI_Comm comm,DA *daref)
1092: {
1093:   int M,N,P,ierr;
1094:   DA  da2;


1099:   if (DAXPeriodic(da->wrap)){
1100:     M = da->refine_x*da->M;
1101:   } else {
1102:     M = da->refine_x*da->M - 1;
1103:   }
1104:   if (DAYPeriodic(da->wrap)){
1105:     N = da->refine_y*da->N;
1106:   } else {
1107:     N = da->refine_y*da->N - 1;
1108:   }
1109:   if (DAZPeriodic(da->wrap)){
1110:     P = da->refine_z*da->P;
1111:   } else {
1112:     P = da->refine_z*da->P - 1;
1113:   }
1114:   if (da->dim == 1) {
1115:     DACreate1d(da->comm,da->wrap,M,da->w,da->s,PETSC_NULL,&da2);
1116:   } else if (da->dim == 2) {
1117:     DACreate2d(da->comm,da->wrap,da->stencil_type,M,N,da->m,da->n,da->w,da->s,PETSC_NULL,PETSC_NULL,&da2);
1118:   } else if (da->dim == 3) {
1119:     DACreate3d(da->comm,da->wrap,da->stencil_type,M,N,P,da->m,da->n,da->p,da->w,da->s,0,0,0,&da2);
1120:   }
1121:   *daref = da2;
1122:   return(0);
1123: }

1125: /*
1126:       M is number of grid points 
1127:       m is number of processors

1129: */
1130: int DASplitComm2d(MPI_Comm comm,int M,int N,int sw,MPI_Comm *outcomm)
1131: {
1132:   int ierr,m,n = 0,csize,size,rank,x = 0,y = 0;

1135:   MPI_Comm_size(comm,&size);
1136:   MPI_Comm_rank(comm,&rank);

1138:   csize = 4*size;
1139:   do {
1140:     if (csize % 4) SETERRQ4(1,"Cannot split communicator of size %d tried %d %d %d",size,csize,x,y);
1141:     csize   = csize/4;
1142: 
1143:     m = (int)(0.5 + sqrt(((double)M)*((double)csize)/((double)N)));
1144:     if (!m) m = 1;
1145:     while (m > 0) {
1146:       n = csize/m;
1147:       if (m*n == csize) break;
1148:       m--;
1149:     }
1150:     if (M > N && m < n) {int _m = m; m = n; n = _m;}

1152:     x = M/m + ((M % m) > ((csize-1) % m));
1153:     y = (N + (csize-1)/m)/n;
1154:   } while ((x < 4 || y < 4) && csize > 1);
1155:   if (size != csize) {
1156:     MPI_Group entire_group,sub_group;
1157:     int       i,*groupies;

1159:     ierr     = MPI_Comm_group(comm,&entire_group);
1160:     PetscMalloc(csize*sizeof(int),&groupies);
1161:     for (i=0; i<csize; i++) {
1162:       groupies[i] = (rank/csize)*csize + i;
1163:     }
1164:     ierr     = MPI_Group_incl(entire_group,csize,groupies,&sub_group);
1165:     ierr     = PetscFree(groupies);
1166:     ierr     = MPI_Comm_create(comm,sub_group,outcomm);
1167:     ierr     = MPI_Group_free(&entire_group);
1168:     ierr     = MPI_Group_free(&sub_group);
1169:     PetscLogInfo(0,"Creating redundant coarse problems of size %dn",csize);
1170:   } else {
1171:     *outcomm = comm;
1172:   }
1173:   return(0);
1174: }